7#include <boost/math/special_functions/pow.hpp>
45template <
int DisplacementDim>
57template <
int DisplacementDim>
61template <
int DisplacementDim>
93 :
value{1 + gamma_p * theta},
104template <
int DisplacementDim>
107 double const sqrtPhi,
double const alpha_p,
double const beta_p,
108 double const delta_p,
double const epsilon_p)
112 4 * boost::math::pow<2>(delta_p) * boost::math::pow<3>(s.
I_1)) /
114 3 * beta_p + 6 * epsilon_p * s.
I_1;
117template <
int DisplacementDim>
122 double const gamma_p,
double const m_p)
125 (s.
D + s.
J_2 * m_p * gamma_p * dtheta_dsigma / one_gt.
value)) /
129template <
int DisplacementDim>
134 double const I_1_squared = boost::math::pow<2>(s.
I_1);
137 return std::sqrt(s.
J_2 * std::pow(1 + mp.
gamma * s.
J_3 /
138 (s.
J_2 * std::sqrt(s.
J_2)),
140 mp.
alpha / 2. * I_1_squared +
141 boost::math::pow<2>(mp.
delta) *
142 boost::math::pow<2>(I_1_squared)) +
146template <
int DisplacementDim>
147std::pair<MathLib::KelvinVector::KelvinVectorType<DisplacementDim>,
152 constexpr int KelvinVectorSize =
162 return {KelvinVector::Zero(), KelvinMatrix::Zero()};
165 auto const& P_dev = Invariants::deviatoric_projection;
168 if (Invariants::determinant(s.
D) == 0)
170 OGS_FATAL(
"Determinant is zero. Matrix is non-invertable.");
174 KelvinVector
const sigma_D_inverse_D = P_dev * sigma_D_inverse;
176 KelvinVector
const dtheta_dsigma =
177 theta * sigma_D_inverse_D - 3. / 2. * theta / s.
J_2 * s.
D;
178 KelvinMatrix
const d2theta_dsigma2 =
180 sigma_D_inverse_D * dtheta_dsigma.transpose() -
181 3. / 2. * theta / s.
J_2 * P_dev -
182 3. / 2. * dtheta_dsigma / s.
J_2 * s.
D.transpose() +
183 3. / 2. * theta / boost::math::pow<2>(s.
J_2) * s.
D * s.
D.transpose();
185 return {dtheta_dsigma, d2theta_dsigma2};
188template <
int DisplacementDim>
196 double const eps_p_V,
197 double const eps_p_V_dot,
198 double const eps_p_eff_dot,
203 static int const KelvinVectorSize =
209 auto const& identity2 = Invariants::identity2;
211 double const theta = s.
J_3 / (s.
J_2 * std::sqrt(s.
J_2));
215 residual.template segment<KelvinVectorSize>(0).noalias() =
216 s.
value / mp.
G - 2 * (eps_D - eps_p_D) -
217 mp.
K / mp.
G * (eps_V - eps_p_V) * identity2;
219 auto const [dtheta_dsigma, d2theta_dsigma2] =
223 double const sqrtPhi = std::sqrt(
225 boost::math::pow<2>(mp.
delta_p) * boost::math::pow<4>(s.
I_1));
227 s, one_gt, sqrtPhi, dtheta_dsigma, mp.
gamma_p, mp.
m_p);
228 KelvinVector
const lambda_flow_D = lambda * flow_D;
230 residual.template segment<KelvinVectorSize>(KelvinVectorSize).noalias() =
231 eps_p_D_dot - lambda_flow_D;
237 residual(2 * KelvinVectorSize, 0) = eps_p_V_dot - lambda * flow_V;
241 residual(2 * KelvinVectorSize + 1) =
243 std::sqrt(2. / 3. * lambda_flow_D.transpose() * lambda_flow_D);
246 residual(2 * KelvinVectorSize + 2) =
yieldFunction(mp, s, k) / mp.
G;
250template <
int DisplacementDim>
257 static int const KelvinVectorSize =
265 auto const& P_dev = Invariants::deviatoric_projection;
266 auto const& identity2 = Invariants::identity2;
268 double const theta = s.
J_3 / (s.
J_2 * std::sqrt(s.
J_2));
271 auto const [dtheta_dsigma, d2theta_dsigma2] =
275 double const sqrtPhi = std::sqrt(
277 boost::math::pow<2>(mp.
delta_p) * boost::math::pow<4>(s.
I_1));
279 s, one_gt, sqrtPhi, dtheta_dsigma, mp.
gamma_p, mp.
m_p);
280 KelvinVector
const lambda_flow_D = lambda * flow_D;
286 jacobian.template block<KelvinVectorSize, KelvinVectorSize>(0, 0)
287 .noalias() = KelvinMatrix::Identity();
291 .template block<KelvinVectorSize, KelvinVectorSize>(0, KelvinVectorSize)
292 .noalias() = 2 * KelvinMatrix::Identity();
295 jacobian.template block<KelvinVectorSize, 1>(0, 2 * KelvinVectorSize)
296 .noalias() = mp.
K / mp.
G * identity2;
304 KelvinVector
const M0 = s.
J_2 / one_gt.
value * dtheta_dsigma;
306 KelvinVector
const dPhi_dsigma =
307 one_gt.
pow_m_p * (s.
D + gm_p * M0) +
309 4 * boost::math::pow<2>(mp.
delta_p) * boost::math::pow<3>(s.
I_1)) *
313 KelvinMatrix
const M1 =
315 (s.
D * dPhi_dsigma.transpose() + gm_p * M0 * dPhi_dsigma.transpose());
317 KelvinMatrix
const M2 =
318 one_gt.
pow_m_p * (P_dev + s.
D * gm_p * M0.transpose());
321 KelvinMatrix
const M3 =
323 ((s.
D + (gm_p - mp.
gamma_p) * M0) * dtheta_dsigma.transpose() +
324 s.
J_2 * d2theta_dsigma2);
327 KelvinMatrix
const dflow_D_dsigma =
328 (-M1 / (4 * boost::math::pow<3>(sqrtPhi)) + (M2 + M3) / (2 * sqrtPhi)) *
331 .template block<KelvinVectorSize, KelvinVectorSize>(KelvinVectorSize, 0)
332 .noalias() = -lambda * dflow_D_dsigma;
336 .template block<KelvinVectorSize, KelvinVectorSize>(KelvinVectorSize,
338 .noalias() = KelvinMatrix::Identity() / dt;
344 .template block<KelvinVectorSize, 1>(KelvinVectorSize,
345 2 * KelvinVectorSize + 2)
346 .noalias() = -flow_D;
351 KelvinVector
const dflow_V_dsigma =
354 boost::math::pow<3>(s.
I_1)) /
355 (4 * boost::math::pow<3>(sqrtPhi)) * dPhi_dsigma +
357 12 * boost::math::pow<2>(mp.
delta_p * s.
I_1) * identity2) /
361 jacobian.template block<1, KelvinVectorSize>(2 * KelvinVectorSize, 0)
362 .noalias() = -lambda * dflow_V_dsigma.transpose();
368 jacobian(2 * KelvinVectorSize, 2 * KelvinVectorSize) = 1. / dt;
376 jacobian(2 * KelvinVectorSize, 2 * KelvinVectorSize + 2) = -flow_V;
380 double const eff_flow =
381 std::sqrt(2. / 3. * lambda_flow_D.transpose() * lambda_flow_D);
386 KelvinVector
const eff_flow23_lambda_flow_D =
387 -2 / 3. / eff_flow * lambda_flow_D;
390 .template block<1, KelvinVectorSize>(2 * KelvinVectorSize + 1, 0)
391 .noalias() = lambda * dflow_D_dsigma * eff_flow23_lambda_flow_D;
393 jacobian(2 * KelvinVectorSize + 1, 2 * KelvinVectorSize + 2) =
394 eff_flow23_lambda_flow_D.transpose() * flow_D;
400 jacobian(2 * KelvinVectorSize + 1, 2 * KelvinVectorSize + 1) = 1. / dt;
404 double const one_gt_pow_m = std::pow(one_gt.
value, mp.
m);
405 double const gm = mp.
gamma * mp.
m;
407 KelvinVector
const dF_dsigma =
409 (one_gt_pow_m * (s.
D + gm * M0) +
411 boost::math::pow<3>(s.
I_1)) *
417 .template block<1, KelvinVectorSize>(2 * KelvinVectorSize + 2, 0)
418 .noalias() = dF_dsigma.transpose() / mp.
G;
422 jacobian(2 * KelvinVectorSize + 2, 2 * KelvinVectorSize + 1) =
431template <
int DisplacementDim>
433 double const K,
double const G)
435 static int const KelvinVectorSize =
439 auto const& P_dev = Invariants::deviatoric_projection;
440 auto const& P_sph = Invariants::spherical_projection;
444 return -2. * I * P_dev - 3. * K / G * I * P_sph;
448 double const hardening_coefficient,
449 double const eps_p_eff)
451 return kappa * (1. + eps_p_eff * hardening_coefficient);
454template <
int DisplacementDim>
456 double const G,
double const K,
462 static int const KelvinVectorSize =
465 auto const& P_dev = Invariants::deviatoric_projection;
468 double const pressure_prev = Invariants::trace(sigma_prev) / (-3. * G);
470 double const e_prev = Invariants::trace(eps_prev);
472 double const pressure = pressure_prev - K / G * (eps_V - e_prev);
475 P_dev * sigma_prev / G;
478 sigma_D_prev + 2 * P_dev * (eps - eps_prev);
479 return sigma_D - pressure * Invariants::identity2;
484template <
typename Res
idualVector,
typename KelvinVector>
485std::tuple<KelvinVector, PlasticStrain<KelvinVector>,
double>
488 static auto const size = KelvinVector::SizeAtCompileTime;
489 return std::forward_as_tuple(
490 solution.template segment<size>(size * 0),
492 solution[size * 2], solution[size * 2 + 1]},
493 solution[size * 2 + 2]);
496template <
int DisplacementDim>
500 std::unique_ptr<DamagePropertiesParameters>&& damage_properties,
503 _mp(std::move(material_properties)),
509template <
int DisplacementDim>
517 material_state_variables)
const
520 &material_state_variables) !=
nullptr);
523 material_state_variables)
526 auto const& identity2 = Invariants::identity2;
527 return (eps - eps_p.D - eps_p.V / 3 * identity2).dot(sigma) / 2;
530template <
int DisplacementDim>
531std::optional<std::tuple<typename SolidEhlers<DisplacementDim>::KelvinVector,
533 DisplacementDim>::MaterialStateVariables>,
540 material_state_variables)
const
542 auto const& eps_m = std::get<MPL::SymmetricTensor<DisplacementDim>>(
544 auto const& eps_m_prev = std::get<MPL::SymmetricTensor<DisplacementDim>>(
546 auto const& sigma_prev = std::get<MPL::SymmetricTensor<DisplacementDim>>(
547 variable_array_prev.
stress);
550 &material_state_variables) !=
nullptr);
554 material_state_variables);
560 double const eps_V = Invariants::trace(eps_m);
562 auto const& P_dev = Invariants::deviatoric_projection;
570 mp.
G, mp.
K, sigma_prev, eps_m, eps_m_prev, eps_V);
580 if ((sigma.squaredNorm() == 0. || dt == 0. ||
584 state.
eps_p.eff)) < 0.))
587 mp.
K - 2. / 3 * mp.
G, mp.
G);
601 Eigen::Matrix<double, JacobianResidualSize, 1>;
615 auto const& eps_p_D =
616 solution.template segment<KelvinVectorSize>(
622 double const eps_p_V_dot = (eps_p_V - state.
eps_p_prev.V) / dt;
625 double const eps_p_eff_dot =
633 solution.template segment<KelvinVectorSize>(
646 auto const update_solution =
649 solution += increment;
651 mp.
G * solution.template segment<KelvinVectorSize>(0)};
655 linear_solver, update_jacobian, update_residual,
658 auto const success_iterations = newton_solver.solve(jacobian);
660 if (!success_iterations)
669 if (*success_iterations == 0)
671 linear_solver.compute(jacobian);
674 std::tie(sigma, state.
eps_p, std::ignore) =
683 Eigen::RowMajor>::Zero();
684 dresidual_deps.template block<KelvinVectorSize, KelvinVectorSize>(0, 0)
697 linear_solver.solve(-dresidual_deps)
698 .template block<KelvinVectorSize, KelvinVectorSize>(0, 0);
707 "Unimplemented tangent type behaviour for the tangent type "
715 return {std::make_tuple(
724template <
int DisplacementDim>
725std::vector<typename MechanicsBase<DisplacementDim>::InternalVariable>
728 return {{
"damage.kappa_d", 1,
731 std::vector<double>& cache) -> std::vector<double>
const&
735 auto const& ehlers_state =
743 state) -> std::span<double>
750 return {&ehlers_state.damage.kappa_d(), 1};
755 std::vector<double>& cache) -> std::vector<double>
const&
759 auto const& ehlers_state =
767 state) -> std::span<double>
774 return {&ehlers_state.damage.value(), 1};
776 {
"eps_p.D", KelvinVector::RowsAtCompileTime,
779 std::vector<double>& cache) -> std::vector<double>
const&
783 auto const& ehlers_state =
786 cache.resize(KelvinVector::RowsAtCompileTime);
788 cache, KelvinVector::RowsAtCompileTime) =
790 ehlers_state.eps_p.D);
795 state) -> std::span<double>
803 ehlers_state.eps_p.D.data(),
804 static_cast<std::size_t
>(KelvinVector::RowsAtCompileTime)};
809 std::vector<double>& cache) -> std::vector<double>
const&
813 auto const& ehlers_state =
817 cache.front() = ehlers_state.
eps_p.V;
821 state) -> std::span<double>
828 return {&ehlers_state.eps_p.V, 1};
833 std::vector<double>& cache) -> std::vector<double>
const&
837 auto const& ehlers_state =
841 cache.front() = ehlers_state.
eps_p.eff;
845 state) -> std::span<double>
852 return {&ehlers_state.eps_p.eff, 1};
868 result(0, 0) = v(0) * v(0);
869 result(0, 1) = result(1, 0) = v(3) * v(3) / 2.;
870 result(0, 2) = result(2, 0) = v(5) * v(5) / 2.;
871 result(0, 3) = result(3, 0) = v(0) * v(3);
872 result(0, 4) = result(4, 0) = v(3) * v(5) / std::sqrt(2.);
873 result(0, 5) = result(5, 0) = v(0) * v(5);
875 result(1, 1) = v(1) * v(1);
876 result(1, 2) = result(2, 1) = v(4) * v(4) / 2.;
877 result(1, 3) = result(3, 1) = v(3) * v(1);
878 result(1, 4) = result(4, 1) = v(1) * v(4);
879 result(1, 5) = result(5, 1) = v(3) * v(4) / std::sqrt(2.);
881 result(2, 2) = v(2) * v(2);
882 result(2, 3) = result(3, 2) = v(5) * v(4) / std::sqrt(2.);
883 result(2, 4) = result(4, 2) = v(4) * v(2);
884 result(2, 5) = result(5, 2) = v(5) * v(2);
886 result(3, 3) = v(0) * v(1) + v(3) * v(3) / 2.;
887 result(3, 4) = result(4, 3) =
888 v(3) * v(4) / 2. + v(5) * v(1) / std::sqrt(2.);
889 result(3, 5) = result(5, 3) =
890 v(0) * v(4) / std::sqrt(2.) + v(3) * v(5) / 2.;
892 result(4, 4) = v(1) * v(2) + v(4) * v(4) / 2.;
893 result(4, 5) = result(5, 4) =
894 v(3) * v(2) / std::sqrt(2.) + v(5) * v(4) / 2.;
896 result(5, 5) = v(0) * v(2) + v(5) * v(5) / 2.;
906 result(0, 0) = v(0) * v(0);
907 result(0, 1) = result(1, 0) = v(3) * v(3) / 2.;
908 result(0, 2) = result(2, 0) = 0;
909 result(0, 3) = result(3, 0) = v(0) * v(3);
911 result(1, 1) = v(1) * v(1);
912 result(1, 2) = result(2, 1) = 0;
913 result(1, 3) = result(3, 1) = v(3) * v(1);
915 result(2, 2) = v(2) * v(2);
916 result(2, 3) = result(3, 2) = 0;
918 result(3, 3) = v(0) * v(1) + v(3) * v(3) / 2.;
std::vector< typename MechanicsBase< DisplacementDim >::InternalVariable > getInternalVariables() const override
MaterialPropertiesParameters _mp
NumLib::NewtonRaphsonSolverParameters const _nonlinear_solver_parameters
static int const KelvinVectorSize
Eigen::Matrix< double, JacobianResidualSize, JacobianResidualSize, Eigen::RowMajor > JacobianMatrix
TangentType const _tangent_type
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
static int const JacobianResidualSize
MathLib::KelvinVector::KelvinVectorType< DisplacementDim > KelvinVector
std::unique_ptr< DamagePropertiesParameters > _damage_properties
MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > KelvinMatrix
Eigen::Matrix< double, JacobianResidualSize, 1 > ResidualVectorType
SolidEhlers(NumLib::NewtonRaphsonSolverParameters nonlinear_solver_parameters, MaterialPropertiesParameters material_properties, std::unique_ptr< DamagePropertiesParameters > &&damage_properties, TangentType tangent_type)
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
KelvinVector mechanical_strain
constexpr auto to_underlying(E e) noexcept
Converts an enumeration to its underlying type.
MathLib::KelvinVector::KelvinMatrixType< 3 > sOdotS< 3 >(MathLib::KelvinVector::KelvinVectorType< 3 > const &v)
MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > calculateDResidualDEps(double const K, double const G)
SolidEhlers< DisplacementDim >::JacobianMatrix calculatePlasticJacobian(double const dt, PhysicalStressWithInvariants< DisplacementDim > const &s, double const lambda, MaterialProperties const &mp)
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)
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)
double calculateIsotropicHardening(double const kappa, double const hardening_coefficient, double const eps_p_eff)
MathLib::KelvinVector::KelvinMatrixType< 2 > sOdotS< 2 >(MathLib::KelvinVector::KelvinVectorType< 2 > const &v)
std::tuple< KelvinVector, PlasticStrain< KelvinVector >, double > splitSolutionVector(ResidualVector const &solution)
std::pair< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > thetaSigmaDerivatives(double theta, PhysicalStressWithInvariants< DisplacementDim > const &s)
double yieldFunction(MaterialProperties const &mp, PhysicalStressWithInvariants< DisplacementDim > const &s, double const k)
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)
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)
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.
double const hardening_coefficient
Holds powers of 1 + gamma_p*theta to base 0, m_p, and m_p-1.
OnePlusGamma_pTheta(double const gamma_p, double const theta, double const m_p)
PhysicalStressWithInvariants(KelvinVector const &stress)
static int const KelvinVectorSize
MathLib::KelvinVector::Invariants< KelvinVectorSize > Invariants
EIGEN_MAKE_ALIGNED_OPERATOR_NEW
MathLib::KelvinVector::KelvinVectorType< DisplacementDim > KelvinVector
void setInitialConditions()
Damage damage
damage part of the state.
PlasticStrain< KelvinVector > eps_p
plastic part of the state.
double getEquivalentPlasticStrain() const override
PlasticStrain< KelvinVector > eps_p_prev
plastic part of the state.
static double FrobeniusNorm(Eigen::Matrix< double, KelvinVectorSize, 1 > const &deviatoric_v)
Get the norm of the deviatoric stress.