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