OGS
Ehlers.cpp
Go to the documentation of this file.
1
10#include "Ehlers.h"
11
12#include <Eigen/LU>
13#include <boost/math/special_functions/pow.hpp>
14
15#include "BaseLib/cpp23.h"
19
20namespace MPL = MaterialPropertyLib;
21
45namespace MaterialLib
46{
47namespace Solids
48{
49namespace Ehlers
50{
51template <int DisplacementDim>
53{
54 return std::sqrt(2.0 / 3.0 * Invariants::FrobeniusNorm(eps_p.D.eval()));
55}
56
63template <int DisplacementDim>
66
67template <int DisplacementDim>
93
96{
97 OnePlusGamma_pTheta(double const gamma_p, double const theta,
98 double const m_p)
99 : value{1 + gamma_p * theta},
100 pow_m_p{std::pow(value, m_p)},
102 {
103 }
104
105 double const value;
106 double const pow_m_p;
107 double const pow_m_p1;
108};
109
110template <int DisplacementDim>
113 double const sqrtPhi, double const alpha_p, double const beta_p,
114 double const delta_p, double const epsilon_p)
115{
116 return 3 *
117 (alpha_p * s.I_1 +
118 4 * boost::math::pow<2>(delta_p) * boost::math::pow<3>(s.I_1)) /
119 (2 * sqrtPhi) +
120 3 * beta_p + 6 * epsilon_p * s.I_1;
121}
122
123template <int DisplacementDim>
126 OnePlusGamma_pTheta const& one_gt, double const sqrtPhi,
127 typename SolidEhlers<DisplacementDim>::KelvinVector const& dtheta_dsigma,
128 double const gamma_p, double const m_p)
129{
130 return (one_gt.pow_m_p *
131 (s.D + s.J_2 * m_p * gamma_p * dtheta_dsigma / one_gt.value)) /
132 (2 * sqrtPhi);
133}
134
135template <int DisplacementDim>
138 double const k)
139{
140 double const I_1_squared = boost::math::pow<2>(s.I_1);
141 assert(s.J_2 != 0);
142
143 return std::sqrt(s.J_2 * std::pow(1 + mp.gamma * s.J_3 /
144 (s.J_2 * std::sqrt(s.J_2)),
145 mp.m) +
146 mp.alpha / 2. * I_1_squared +
147 boost::math::pow<2>(mp.delta) *
148 boost::math::pow<2>(I_1_squared)) +
149 mp.beta * s.I_1 + mp.epsilon * I_1_squared - k;
150}
151
152template <int DisplacementDim>
153std::pair<MathLib::KelvinVector::KelvinVectorType<DisplacementDim>,
157{
158 constexpr int KelvinVectorSize =
161 using KelvinVector =
163 using KelvinMatrix =
165
166 if (theta == 0)
167 {
168 return {KelvinVector::Zero(), KelvinMatrix::Zero()};
169 }
170
171 auto const& P_dev = Invariants::deviatoric_projection;
172
173 // inverse of deviatoric stress tensor
174 if (Invariants::determinant(s.D) == 0)
175 {
176 OGS_FATAL("Determinant is zero. Matrix is non-invertable.");
177 }
178 // inverse of sigma_D
179 KelvinVector const sigma_D_inverse = MathLib::KelvinVector::inverse(s.D);
180 KelvinVector const sigma_D_inverse_D = P_dev * sigma_D_inverse;
181
182 KelvinVector const dtheta_dsigma =
183 theta * sigma_D_inverse_D - 3. / 2. * theta / s.J_2 * s.D;
184 KelvinMatrix const d2theta_dsigma2 =
185 theta * P_dev * sOdotS<DisplacementDim>(sigma_D_inverse) * P_dev +
186 sigma_D_inverse_D * dtheta_dsigma.transpose() -
187 3. / 2. * theta / s.J_2 * P_dev -
188 3. / 2. * dtheta_dsigma / s.J_2 * s.D.transpose() +
189 3. / 2. * theta / boost::math::pow<2>(s.J_2) * s.D * s.D.transpose();
190
191 return {dtheta_dsigma, d2theta_dsigma2};
192}
193
194template <int DisplacementDim>
198 double const eps_V,
202 double const eps_p_V,
203 double const eps_p_V_dot,
204 double const eps_p_eff_dot,
205 double const lambda,
206 double const k,
207 MaterialProperties const& mp)
208{
209 static int const KelvinVectorSize =
212 using KelvinVector =
214
215 auto const& identity2 = Invariants::identity2;
216
217 double const theta = s.J_3 / (s.J_2 * std::sqrt(s.J_2));
218
220 // calculate stress residual
221 residual.template segment<KelvinVectorSize>(0).noalias() =
222 s.value / mp.G - 2 * (eps_D - eps_p_D) -
223 mp.K / mp.G * (eps_V - eps_p_V) * identity2;
224
225 auto const [dtheta_dsigma, d2theta_dsigma2] =
226 thetaSigmaDerivatives(theta, s);
227
228 OnePlusGamma_pTheta const one_gt{mp.gamma_p, theta, mp.m_p};
229 double const sqrtPhi = std::sqrt(
230 s.J_2 * one_gt.pow_m_p + mp.alpha_p / 2. * boost::math::pow<2>(s.I_1) +
231 boost::math::pow<2>(mp.delta_p) * boost::math::pow<4>(s.I_1));
232 KelvinVector const flow_D = plasticFlowDeviatoricPart(
233 s, one_gt, sqrtPhi, dtheta_dsigma, mp.gamma_p, mp.m_p);
234 KelvinVector const lambda_flow_D = lambda * flow_D;
235
236 residual.template segment<KelvinVectorSize>(KelvinVectorSize).noalias() =
237 eps_p_D_dot - lambda_flow_D;
238
239 // plastic volume strain
240 {
241 double const flow_V = plasticFlowVolumetricPart<DisplacementDim>(
242 s, sqrtPhi, mp.alpha_p, mp.beta_p, mp.delta_p, mp.epsilon_p);
243 residual(2 * KelvinVectorSize, 0) = eps_p_V_dot - lambda * flow_V;
244 }
245
246 // evolution of plastic equivalent strain
247 residual(2 * KelvinVectorSize + 1) =
248 eps_p_eff_dot -
249 std::sqrt(2. / 3. * lambda_flow_D.transpose() * lambda_flow_D);
250
251 // yield function (for plastic multiplier)
252 residual(2 * KelvinVectorSize + 2) = yieldFunction(mp, s, k) / mp.G;
253 return residual;
254}
255
256template <int DisplacementDim>
258 double const dt,
260 double const lambda,
261 MaterialProperties const& mp)
262{
263 static int const KelvinVectorSize =
266 using KelvinVector =
268 using KelvinMatrix =
270
271 auto const& P_dev = Invariants::deviatoric_projection;
272 auto const& identity2 = Invariants::identity2;
273
274 double const theta = s.J_3 / (s.J_2 * std::sqrt(s.J_2));
275 OnePlusGamma_pTheta const one_gt{mp.gamma_p, theta, mp.m_p};
276
277 auto const [dtheta_dsigma, d2theta_dsigma2] =
278 thetaSigmaDerivatives(theta, s);
279
280 // deviatoric flow
281 double const sqrtPhi = std::sqrt(
282 s.J_2 * one_gt.pow_m_p + mp.alpha_p / 2. * boost::math::pow<2>(s.I_1) +
283 boost::math::pow<2>(mp.delta_p) * boost::math::pow<4>(s.I_1));
284 KelvinVector const flow_D = plasticFlowDeviatoricPart(
285 s, one_gt, sqrtPhi, dtheta_dsigma, mp.gamma_p, mp.m_p);
286 KelvinVector const lambda_flow_D = lambda * flow_D;
287
290
291 // G_11
292 jacobian.template block<KelvinVectorSize, KelvinVectorSize>(0, 0)
293 .noalias() = KelvinMatrix::Identity();
294
295 // G_12
296 jacobian
297 .template block<KelvinVectorSize, KelvinVectorSize>(0, KelvinVectorSize)
298 .noalias() = 2 * KelvinMatrix::Identity();
299
300 // G_13
301 jacobian.template block<KelvinVectorSize, 1>(0, 2 * KelvinVectorSize)
302 .noalias() = mp.K / mp.G * identity2;
303
304 // G_14 and G_15 are zero
305
306 // G_21 -- derivative of deviatoric flow
307
308 double const gm_p = mp.gamma_p * mp.m_p;
309 // intermediate variable for derivative of deviatoric flow
310 KelvinVector const M0 = s.J_2 / one_gt.value * dtheta_dsigma;
311 // derivative of Phi w.r.t. sigma
312 KelvinVector const dPhi_dsigma =
313 one_gt.pow_m_p * (s.D + gm_p * M0) +
314 (mp.alpha_p * s.I_1 +
315 4 * boost::math::pow<2>(mp.delta_p) * boost::math::pow<3>(s.I_1)) *
316 identity2;
317
318 // intermediate variable for derivative of deviatoric flow
319 KelvinMatrix const M1 =
320 one_gt.pow_m_p *
321 (s.D * dPhi_dsigma.transpose() + gm_p * M0 * dPhi_dsigma.transpose());
322 // intermediate variable for derivative of deviatoric flow
323 KelvinMatrix const M2 =
324 one_gt.pow_m_p * (P_dev + s.D * gm_p * M0.transpose());
325
326 // intermediate variable for derivative of deviatoric flow
327 KelvinMatrix const M3 =
328 gm_p * one_gt.pow_m_p1 *
329 ((s.D + (gm_p - mp.gamma_p) * M0) * dtheta_dsigma.transpose() +
330 s.J_2 * d2theta_dsigma2);
331
332 // derivative of flow_D w.r.t. sigma
333 KelvinMatrix const dflow_D_dsigma =
334 (-M1 / (4 * boost::math::pow<3>(sqrtPhi)) + (M2 + M3) / (2 * sqrtPhi)) *
335 mp.G;
336 jacobian
337 .template block<KelvinVectorSize, KelvinVectorSize>(KelvinVectorSize, 0)
338 .noalias() = -lambda * dflow_D_dsigma;
339
340 // G_22
341 jacobian
342 .template block<KelvinVectorSize, KelvinVectorSize>(KelvinVectorSize,
343 KelvinVectorSize)
344 .noalias() = KelvinMatrix::Identity() / dt;
345
346 // G_23 and G_24 are zero
347
348 // G_25
349 jacobian
350 .template block<KelvinVectorSize, 1>(KelvinVectorSize,
351 2 * KelvinVectorSize + 2)
352 .noalias() = -flow_D;
353
354 // G_31
355 {
356 // derivative of flow_V w.r.t. sigma
357 KelvinVector const dflow_V_dsigma =
358 3 * mp.G *
359 (-(mp.alpha_p * s.I_1 + 4 * boost::math::pow<2>(mp.delta_p) *
360 boost::math::pow<3>(s.I_1)) /
361 (4 * boost::math::pow<3>(sqrtPhi)) * dPhi_dsigma +
362 (mp.alpha_p * identity2 +
363 12 * boost::math::pow<2>(mp.delta_p * s.I_1) * identity2) /
364 (2 * sqrtPhi) +
365 2 * mp.epsilon_p * identity2);
366
367 jacobian.template block<1, KelvinVectorSize>(2 * KelvinVectorSize, 0)
368 .noalias() = -lambda * dflow_V_dsigma.transpose();
369 }
370
371 // G_32 is zero
372
373 // G_33
374 jacobian(2 * KelvinVectorSize, 2 * KelvinVectorSize) = 1. / dt;
375
376 // G_34 is zero
377
378 // G_35
379 {
380 double const flow_V = plasticFlowVolumetricPart<DisplacementDim>(
381 s, sqrtPhi, mp.alpha_p, mp.beta_p, mp.delta_p, mp.epsilon_p);
382 jacobian(2 * KelvinVectorSize, 2 * KelvinVectorSize + 2) = -flow_V;
383 }
384
385 // increment of effectiv plastic strain
386 double const eff_flow =
387 std::sqrt(2. / 3. * lambda_flow_D.transpose() * lambda_flow_D);
388
389 if (eff_flow > 0)
390 {
391 // intermediate variable for derivative of plastic jacobian
392 KelvinVector const eff_flow23_lambda_flow_D =
393 -2 / 3. / eff_flow * lambda_flow_D;
394 // G_41
395 jacobian
396 .template block<1, KelvinVectorSize>(2 * KelvinVectorSize + 1, 0)
397 .noalias() = lambda * dflow_D_dsigma * eff_flow23_lambda_flow_D;
398 // G_45
399 jacobian(2 * KelvinVectorSize + 1, 2 * KelvinVectorSize + 2) =
400 eff_flow23_lambda_flow_D.transpose() * flow_D;
401 }
402
403 // G_42 and G_43 are zero
404
405 // G_44
406 jacobian(2 * KelvinVectorSize + 1, 2 * KelvinVectorSize + 1) = 1. / dt;
407
408 // G_51
409 {
410 double const one_gt_pow_m = std::pow(one_gt.value, mp.m);
411 double const gm = mp.gamma * mp.m;
412 // derivative of yield function w.r.t. sigma
413 KelvinVector const dF_dsigma =
414 mp.G *
415 (one_gt_pow_m * (s.D + gm * M0) +
416 (mp.alpha * s.I_1 + 4 * boost::math::pow<2>(mp.delta) *
417 boost::math::pow<3>(s.I_1)) *
418 identity2) /
419 (2. * sqrtPhi) +
420 mp.G * (mp.beta + 2 * mp.epsilon_p * s.I_1) * identity2;
421
422 jacobian
423 .template block<1, KelvinVectorSize>(2 * KelvinVectorSize + 2, 0)
424 .noalias() = dF_dsigma.transpose() / mp.G;
425 }
426
427 // G_54
428 jacobian(2 * KelvinVectorSize + 2, 2 * KelvinVectorSize + 1) =
429 -mp.kappa * mp.hardening_coefficient / mp.G;
430
431 // G_52, G_53, G_55 are zero
432 return jacobian;
433}
434
437template <int DisplacementDim>
439 double const K, double const G)
440{
441 static int const KelvinVectorSize =
444
445 auto const& P_dev = Invariants::deviatoric_projection;
446 auto const& P_sph = Invariants::spherical_projection;
447 auto const& I =
449
450 return -2. * I * P_dev - 3. * K / G * I * P_sph;
451}
452
453inline double calculateIsotropicHardening(double const kappa,
454 double const hardening_coefficient,
455 double const eps_p_eff)
456{
457 return kappa * (1. + eps_p_eff * hardening_coefficient);
458}
459
460template <int DisplacementDim>
462 double const G, double const K,
463 typename SolidEhlers<DisplacementDim>::KelvinVector const& sigma_prev,
465 typename SolidEhlers<DisplacementDim>::KelvinVector const& eps_prev,
466 double const eps_V)
467{
468 static int const KelvinVectorSize =
471 auto const& P_dev = Invariants::deviatoric_projection;
472
473 // dimensionless initial hydrostatic pressure
474 double const pressure_prev = Invariants::trace(sigma_prev) / (-3. * G);
475 // initial strain invariant
476 double const e_prev = Invariants::trace(eps_prev);
477 // dimensioness hydrostatic stress increment
478 double const pressure = pressure_prev - K / G * (eps_V - e_prev);
479 // dimensionless deviatoric initial stress
480 typename SolidEhlers<DisplacementDim>::KelvinVector const sigma_D_prev =
481 P_dev * sigma_prev / G;
482 // dimensionless deviatoric stress
483 typename SolidEhlers<DisplacementDim>::KelvinVector const sigma_D =
484 sigma_D_prev + 2 * P_dev * (eps - eps_prev);
485 return sigma_D - pressure * Invariants::identity2;
486}
487
490template <typename ResidualVector, typename KelvinVector>
491std::tuple<KelvinVector, PlasticStrain<KelvinVector>, double>
492splitSolutionVector(ResidualVector const& solution)
493{
494 static auto const size = KelvinVector::SizeAtCompileTime;
495 return std::forward_as_tuple(
496 solution.template segment<size>(size * 0),
497 PlasticStrain<KelvinVector>{solution.template segment<size>(size * 1),
498 solution[size * 2], solution[size * 2 + 1]},
499 solution[size * 2 + 2]);
500}
501
502template <int DisplacementDim>
504 NumLib::NewtonRaphsonSolverParameters nonlinear_solver_parameters,
505 MaterialPropertiesParameters material_properties,
506 std::unique_ptr<DamagePropertiesParameters>&& damage_properties,
507 TangentType tangent_type)
508 : _nonlinear_solver_parameters(std::move(nonlinear_solver_parameters)),
509 _mp(std::move(material_properties)),
510 _damage_properties(std::move(damage_properties)),
511 _tangent_type(tangent_type)
512{
513}
514
515template <int DisplacementDim>
517 double const /*t*/,
519 double const /*dt*/,
520 KelvinVector const& eps,
521 KelvinVector const& sigma,
523 material_state_variables) const
524{
525 assert(dynamic_cast<StateVariables<DisplacementDim> const*>(
526 &material_state_variables) != nullptr);
527
528 auto const& eps_p = static_cast<StateVariables<DisplacementDim> const&>(
529 material_state_variables)
530 .eps_p;
532 auto const& identity2 = Invariants::identity2;
533 return (eps - eps_p.D - eps_p.V / 3 * identity2).dot(sigma) / 2;
534}
535
536template <int DisplacementDim>
537std::optional<std::tuple<typename SolidEhlers<DisplacementDim>::KelvinVector,
538 std::unique_ptr<typename MechanicsBase<
539 DisplacementDim>::MaterialStateVariables>,
542 MaterialPropertyLib::VariableArray const& variable_array_prev,
543 MaterialPropertyLib::VariableArray const& variable_array, double const t,
544 ParameterLib::SpatialPosition const& x, double const dt,
546 material_state_variables) const
547{
548 auto const& eps_m = std::get<MPL::SymmetricTensor<DisplacementDim>>(
549 variable_array.mechanical_strain);
550 auto const& eps_m_prev = std::get<MPL::SymmetricTensor<DisplacementDim>>(
551 variable_array_prev.mechanical_strain);
552 auto const& sigma_prev = std::get<MPL::SymmetricTensor<DisplacementDim>>(
553 variable_array_prev.stress);
554
555 assert(dynamic_cast<StateVariables<DisplacementDim> const*>(
556 &material_state_variables) != nullptr);
557
559 static_cast<StateVariables<DisplacementDim> const&>(
560 material_state_variables);
561 state.setInitialConditions();
562
564
565 // volumetric strain
566 double const eps_V = Invariants::trace(eps_m);
567
568 auto const& P_dev = Invariants::deviatoric_projection;
569 // deviatoric strain
570 KelvinVector const eps_D = P_dev * eps_m;
571
572 // do the evaluation once per function call.
573 MaterialProperties const mp(t, x, _mp);
574
575 KelvinVector sigma = predict_sigma<DisplacementDim>(
576 mp.G, mp.K, sigma_prev, eps_m, eps_m_prev, eps_V);
577
578 KelvinMatrix tangentStiffness;
579
581 // Quit early if sigma is zero (nothing to do), or dt is zero, or if we are
582 // still in elastic zone.
583 // dt can be zero if we are in the initialization phase and the tangent
584 // stiffness will no be necessary. Anyway the newton loop below would not
585 // work because of division by zero.
586 if ((sigma.squaredNorm() == 0. || dt == 0. ||
588 mp, s,
590 state.eps_p.eff)) < 0.))
591 {
592 tangentStiffness = elasticTangentStiffness<DisplacementDim>(
593 mp.K - 2. / 3 * mp.G, mp.G);
594 }
595 else
596 {
597 // Linear solver for the newton loop is required after the loop with the
598 // same matrix. This saves one decomposition.
599 Eigen::FullPivLU<Eigen::Matrix<double, JacobianResidualSize,
600 JacobianResidualSize, Eigen::RowMajor>>
601 linear_solver;
602
603 {
604 using KelvinVector =
606 using ResidualVectorType =
607 Eigen::Matrix<double, JacobianResidualSize, 1>;
608 using JacobianMatrix =
609 Eigen::Matrix<double, JacobianResidualSize,
610 JacobianResidualSize, Eigen::RowMajor>;
611
612 JacobianMatrix jacobian;
613
614 // Agglomerated solution vector construction. It is later split
615 // into individual parts by splitSolutionVector().
616 ResidualVectorType solution;
617 solution << sigma, state.eps_p.D, state.eps_p.V, state.eps_p.eff, 0;
618
619 auto const update_residual = [&](ResidualVectorType& residual)
620 {
621 auto const& eps_p_D =
622 solution.template segment<KelvinVectorSize>(
623 KelvinVectorSize);
624 KelvinVector const eps_p_D_dot =
625 (eps_p_D - state.eps_p_prev.D) / dt;
626
627 double const& eps_p_V = solution[KelvinVectorSize * 2];
628 double const eps_p_V_dot = (eps_p_V - state.eps_p_prev.V) / dt;
629
630 double const& eps_p_eff = solution[KelvinVectorSize * 2 + 1];
631 double const eps_p_eff_dot =
632 (eps_p_eff - state.eps_p_prev.eff) / dt;
633
634 double const k_hardening = calculateIsotropicHardening(
636 solution[KelvinVectorSize * 2 + 1]);
637 residual = calculatePlasticResidual<DisplacementDim>(
638 eps_D, eps_V, s,
639 solution.template segment<KelvinVectorSize>(
640 KelvinVectorSize),
641 eps_p_D_dot, solution[KelvinVectorSize * 2], eps_p_V_dot,
642 eps_p_eff_dot, solution[KelvinVectorSize * 2 + 2],
643 k_hardening, mp);
644 };
645
646 auto const update_jacobian = [&](JacobianMatrix& jacobian)
647 {
648 jacobian = calculatePlasticJacobian<DisplacementDim>(
649 dt, s, solution[KelvinVectorSize * 2 + 2], mp);
650 };
651
652 auto const update_solution =
653 [&](ResidualVectorType const& increment)
654 {
655 solution += increment;
657 mp.G * solution.template segment<KelvinVectorSize>(0)};
658 };
659
660 auto newton_solver = NumLib::NewtonRaphson<
661 decltype(linear_solver), JacobianMatrix,
662 decltype(update_jacobian), ResidualVectorType,
663 decltype(update_residual), decltype(update_solution)>(
664 linear_solver, update_jacobian, update_residual,
665 update_solution, _nonlinear_solver_parameters);
666
667 auto const success_iterations = newton_solver.solve(jacobian);
668
669 if (!success_iterations)
670 {
671 return {};
672 }
673
674 // If the Newton loop didn't run, the linear solver will not be
675 // initialized.
676 // This happens usually for the first iteration of the first
677 // timestep.
678 if (*success_iterations == 0)
679 {
680 linear_solver.compute(jacobian);
681 }
682
683 std::tie(sigma, state.eps_p, std::ignore) =
684 splitSolutionVector<ResidualVectorType, KelvinVector>(solution);
685 }
686
687 // Calculate residual derivative w.r.t. strain
688 Eigen::Matrix<double, JacobianResidualSize, KelvinVectorSize,
689 Eigen::RowMajor>
690 dresidual_deps =
691 Eigen::Matrix<double, JacobianResidualSize, KelvinVectorSize,
692 Eigen::RowMajor>::Zero();
693 dresidual_deps.template block<KelvinVectorSize, KelvinVectorSize>(0, 0)
694 .noalias() = calculateDResidualDEps<DisplacementDim>(mp.K, mp.G);
695
696 if (_tangent_type == TangentType::Elastic)
697 {
698 tangentStiffness =
699 elasticTangentStiffness<DisplacementDim>(mp.K, mp.G);
700 }
701 else if (_tangent_type == TangentType::Plastic ||
702 _tangent_type == TangentType::PlasticDamageSecant)
703 {
704 tangentStiffness =
705 mp.G *
706 linear_solver.solve(-dresidual_deps)
707 .template block<KelvinVectorSize, KelvinVectorSize>(0, 0);
708 if (_tangent_type == TangentType::PlasticDamageSecant)
709 {
710 tangentStiffness *= 1 - state.damage.value();
711 }
712 }
713 else
714 {
715 OGS_FATAL(
716 "Unimplemented tangent type behaviour for the tangent type "
717 "'{}'.",
718 BaseLib::to_underlying(_tangent_type));
719 }
720 }
721
722 KelvinVector sigma_final = mp.G * sigma;
723
724 return {std::make_tuple(
725 sigma_final,
726 std::unique_ptr<
729 static_cast<StateVariables<DisplacementDim> const&>(state)}},
730 tangentStiffness)};
731}
732
733template <int DisplacementDim>
734std::vector<typename MechanicsBase<DisplacementDim>::InternalVariable>
736{
737 return {{"damage.kappa_d", 1,
738 [](typename MechanicsBase<
739 DisplacementDim>::MaterialStateVariables const& state,
740 std::vector<double>& cache) -> std::vector<double> const&
741 {
742 assert(dynamic_cast<StateVariables<DisplacementDim> const*>(
743 &state) != nullptr);
744 auto const& ehlers_state =
745 static_cast<StateVariables<DisplacementDim> const&>(state);
746
747 cache.resize(1);
748 cache.front() = ehlers_state.damage.kappa_d();
749 return cache;
750 },
752 state) -> std::span<double>
753 {
754 assert(dynamic_cast<StateVariables<DisplacementDim> const*>(
755 &state) != nullptr);
756 auto& ehlers_state =
757 static_cast<StateVariables<DisplacementDim>&>(state);
758
759 return {&ehlers_state.damage.kappa_d(), 1};
760 }},
761 {"damage.value", 1,
762 [](typename MechanicsBase<
763 DisplacementDim>::MaterialStateVariables const& state,
764 std::vector<double>& cache) -> std::vector<double> const&
765 {
766 assert(dynamic_cast<StateVariables<DisplacementDim> const*>(
767 &state) != nullptr);
768 auto const& ehlers_state =
769 static_cast<StateVariables<DisplacementDim> const&>(state);
770
771 cache.resize(1);
772 cache.front() = ehlers_state.damage.value();
773 return cache;
774 },
776 state) -> std::span<double>
777 {
778 assert(dynamic_cast<StateVariables<DisplacementDim> const*>(
779 &state) != nullptr);
780 auto& ehlers_state =
781 static_cast<StateVariables<DisplacementDim>&>(state);
782
783 return {&ehlers_state.damage.value(), 1};
784 }},
785 {"eps_p.D", KelvinVector::RowsAtCompileTime,
786 [](typename MechanicsBase<
787 DisplacementDim>::MaterialStateVariables const& state,
788 std::vector<double>& cache) -> std::vector<double> const&
789 {
790 assert(dynamic_cast<StateVariables<DisplacementDim> const*>(
791 &state) != nullptr);
792 auto const& ehlers_state =
793 static_cast<StateVariables<DisplacementDim> const&>(state);
794
795 cache.resize(KelvinVector::RowsAtCompileTime);
796 MathLib::toVector<KelvinVector>(
797 cache, KelvinVector::RowsAtCompileTime) =
799 ehlers_state.eps_p.D);
800
801 return cache;
802 },
804 state) -> std::span<double>
805 {
806 assert(dynamic_cast<StateVariables<DisplacementDim> const*>(
807 &state) != nullptr);
808 auto& ehlers_state =
809 static_cast<StateVariables<DisplacementDim>&>(state);
810
811 return {
812 ehlers_state.eps_p.D.data(),
813 static_cast<std::size_t>(KelvinVector::RowsAtCompileTime)};
814 }},
815 {"eps_p.V", 1,
816 [](typename MechanicsBase<
817 DisplacementDim>::MaterialStateVariables const& state,
818 std::vector<double>& cache) -> std::vector<double> const&
819 {
820 assert(dynamic_cast<StateVariables<DisplacementDim> const*>(
821 &state) != nullptr);
822 auto const& ehlers_state =
823 static_cast<StateVariables<DisplacementDim> const&>(state);
824
825 cache.resize(1);
826 cache.front() = ehlers_state.eps_p.V;
827 return cache;
828 },
830 state) -> std::span<double>
831 {
832 assert(dynamic_cast<StateVariables<DisplacementDim> const*>(
833 &state) != nullptr);
834 auto& ehlers_state =
835 static_cast<StateVariables<DisplacementDim>&>(state);
836
837 return {&ehlers_state.eps_p.V, 1};
838 }},
839 {"eps_p.eff", 1,
840 [](typename MechanicsBase<
841 DisplacementDim>::MaterialStateVariables const& state,
842 std::vector<double>& cache) -> std::vector<double> const&
843 {
844 assert(dynamic_cast<StateVariables<DisplacementDim> const*>(
845 &state) != nullptr);
846 auto const& ehlers_state =
847 static_cast<StateVariables<DisplacementDim> const&>(state);
848
849 cache.resize(1);
850 cache.front() = ehlers_state.eps_p.eff;
851 return cache;
852 },
854 state) -> std::span<double>
855 {
856 assert(dynamic_cast<StateVariables<DisplacementDim> const*>(
857 &state) != nullptr);
858 auto& ehlers_state =
859 static_cast<StateVariables<DisplacementDim>&>(state);
860
861 return {&ehlers_state.eps_p.eff, 1};
862 }}};
863}
864
865template class SolidEhlers<2>;
866template class SolidEhlers<3>;
867
868template struct StateVariables<2>;
869template struct StateVariables<3>;
870
871template <>
874{
876
877 result(0, 0) = v(0) * v(0);
878 result(0, 1) = result(1, 0) = v(3) * v(3) / 2.;
879 result(0, 2) = result(2, 0) = v(5) * v(5) / 2.;
880 result(0, 3) = result(3, 0) = v(0) * v(3);
881 result(0, 4) = result(4, 0) = v(3) * v(5) / std::sqrt(2.);
882 result(0, 5) = result(5, 0) = v(0) * v(5);
883
884 result(1, 1) = v(1) * v(1);
885 result(1, 2) = result(2, 1) = v(4) * v(4) / 2.;
886 result(1, 3) = result(3, 1) = v(3) * v(1);
887 result(1, 4) = result(4, 1) = v(1) * v(4);
888 result(1, 5) = result(5, 1) = v(3) * v(4) / std::sqrt(2.);
889
890 result(2, 2) = v(2) * v(2);
891 result(2, 3) = result(3, 2) = v(5) * v(4) / std::sqrt(2.);
892 result(2, 4) = result(4, 2) = v(4) * v(2);
893 result(2, 5) = result(5, 2) = v(5) * v(2);
894
895 result(3, 3) = v(0) * v(1) + v(3) * v(3) / 2.;
896 result(3, 4) = result(4, 3) =
897 v(3) * v(4) / 2. + v(5) * v(1) / std::sqrt(2.);
898 result(3, 5) = result(5, 3) =
899 v(0) * v(4) / std::sqrt(2.) + v(3) * v(5) / 2.;
900
901 result(4, 4) = v(1) * v(2) + v(4) * v(4) / 2.;
902 result(4, 5) = result(5, 4) =
903 v(3) * v(2) / std::sqrt(2.) + v(5) * v(4) / 2.;
904
905 result(5, 5) = v(0) * v(2) + v(5) * v(5) / 2.;
906 return result;
907}
908
909template <>
912{
914
915 result(0, 0) = v(0) * v(0);
916 result(0, 1) = result(1, 0) = v(3) * v(3) / 2.;
917 result(0, 2) = result(2, 0) = 0;
918 result(0, 3) = result(3, 0) = v(0) * v(3);
919
920 result(1, 1) = v(1) * v(1);
921 result(1, 2) = result(2, 1) = 0;
922 result(1, 3) = result(3, 1) = v(3) * v(1);
923
924 result(2, 2) = v(2) * v(2);
925 result(2, 3) = result(3, 2) = 0;
926
927 result(3, 3) = v(0) * v(1) + v(3) * v(3) / 2.;
928
929 return result;
930}
931
932} // namespace Ehlers
933} // namespace Solids
934} // namespace MaterialLib
#define OGS_FATAL(...)
Definition Error.h:26
std::vector< typename MechanicsBase< DisplacementDim >::InternalVariable > getInternalVariables() const override
Definition Ehlers.cpp:735
Eigen::Matrix< double, JacobianResidualSize, JacobianResidualSize, Eigen::RowMajor > JacobianMatrix
Definition Ehlers.h:283
std::optional< std::tuple< KelvinVector, std::unique_ptr< typename MechanicsBase< DisplacementDim >::MaterialStateVariables >, KelvinMatrix > > integrateStress(MaterialPropertyLib::VariableArray const &variable_array_prev, MaterialPropertyLib::VariableArray const &variable_array, double const t, ParameterLib::SpatialPosition const &x, double const dt, typename MechanicsBase< DisplacementDim >::MaterialStateVariables const &material_state_variables) const override
Definition Ehlers.cpp:541
MathLib::KelvinVector::KelvinVectorType< DisplacementDim > KelvinVector
Definition Ehlers.h:295
MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > KelvinMatrix
Definition Ehlers.h:297
Eigen::Matrix< double, JacobianResidualSize, 1 > ResidualVectorType
Definition Ehlers.h:282
SolidEhlers(NumLib::NewtonRaphsonSolverParameters nonlinear_solver_parameters, MaterialPropertiesParameters material_properties, std::unique_ptr< DamagePropertiesParameters > &&damage_properties, TangentType tangent_type)
Definition Ehlers.cpp:503
double computeFreeEnergyDensity(double const, ParameterLib::SpatialPosition const &, double const, KelvinVector const &eps, KelvinVector const &sigma, typename MechanicsBase< DisplacementDim >::MaterialStateVariables const &material_state_variables) const override
Definition Ehlers.cpp:516
std::variant< std::monostate, Eigen::Matrix< double, 4, 1 >, Eigen::Matrix< double, 6, 1 > > mechanical_strain
std::variant< std::monostate, Eigen::Matrix< double, 4, 1 >, Eigen::Matrix< double, 6, 1 > > stress
std::optional< int > solve(JacobianMatrix &jacobian) const
constexpr auto to_underlying(E e) noexcept
Converts an enumeration to its underlying type.
Definition cpp23.h:29
MathLib::KelvinVector::KelvinMatrixType< 3 > sOdotS< 3 >(MathLib::KelvinVector::KelvinVectorType< 3 > const &v)
Definition Ehlers.cpp:872
MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > calculateDResidualDEps(double const K, double const G)
Definition Ehlers.cpp:438
SolidEhlers< DisplacementDim >::JacobianMatrix calculatePlasticJacobian(double const dt, PhysicalStressWithInvariants< DisplacementDim > const &s, double const lambda, MaterialProperties const &mp)
Definition Ehlers.cpp:257
SolidEhlers< DisplacementDim >::KelvinVector plasticFlowDeviatoricPart(PhysicalStressWithInvariants< DisplacementDim > const &s, OnePlusGamma_pTheta const &one_gt, double const sqrtPhi, typename SolidEhlers< DisplacementDim >::KelvinVector const &dtheta_dsigma, double const gamma_p, double const m_p)
Definition Ehlers.cpp:124
MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > sOdotS(MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &v)
double plasticFlowVolumetricPart(PhysicalStressWithInvariants< DisplacementDim > const &s, double const sqrtPhi, double const alpha_p, double const beta_p, double const delta_p, double const epsilon_p)
Definition Ehlers.cpp:111
double calculateIsotropicHardening(double const kappa, double const hardening_coefficient, double const eps_p_eff)
Definition Ehlers.cpp:453
MathLib::KelvinVector::KelvinMatrixType< 2 > sOdotS< 2 >(MathLib::KelvinVector::KelvinVectorType< 2 > const &v)
Definition Ehlers.cpp:910
std::tuple< KelvinVector, PlasticStrain< KelvinVector >, double > splitSolutionVector(ResidualVector const &solution)
Definition Ehlers.cpp:492
std::pair< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > thetaSigmaDerivatives(double theta, PhysicalStressWithInvariants< DisplacementDim > const &s)
Definition Ehlers.cpp:155
double yieldFunction(MaterialProperties const &mp, PhysicalStressWithInvariants< DisplacementDim > const &s, double const k)
Definition Ehlers.cpp:136
SolidEhlers< DisplacementDim >::KelvinVector predict_sigma(double const G, double const K, typename SolidEhlers< DisplacementDim >::KelvinVector const &sigma_prev, typename SolidEhlers< DisplacementDim >::KelvinVector const &eps, typename SolidEhlers< DisplacementDim >::KelvinVector const &eps_prev, double const eps_V)
Definition Ehlers.cpp:461
SolidEhlers< DisplacementDim >::ResidualVectorType calculatePlasticResidual(MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps_D, double const eps_V, PhysicalStressWithInvariants< DisplacementDim > const &s, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps_p_D, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps_p_D_dot, double const eps_p_V, double const eps_p_V_dot, double const eps_p_eff_dot, double const lambda, double const k, MaterialProperties const &mp)
Definition Ehlers.cpp:196
Eigen::Matrix< double, 4, 1 > kelvinVectorToSymmetricTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor > KelvinMatrixType
Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > inverse(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
Holds powers of 1 + gamma_p*theta to base 0, m_p, and m_p-1.
Definition Ehlers.cpp:96
OnePlusGamma_pTheta(double const gamma_p, double const theta, double const m_p)
Definition Ehlers.cpp:97
PhysicalStressWithInvariants(KelvinVector const &stress)
Definition Ehlers.cpp:76
MathLib::KelvinVector::KelvinVectorType< DisplacementDim > KelvinVector
Definition Ehlers.cpp:73
Damage damage
damage part of the state.
Definition Ehlers.h:245
PlasticStrain< KelvinVector > eps_p
plastic part of the state.
Definition Ehlers.h:244
double getEquivalentPlasticStrain() const override
Definition Ehlers.cpp:52
PlasticStrain< KelvinVector > eps_p_prev
plastic part of the state.
Definition Ehlers.h:248