7#include <boost/math/special_functions/pow.hpp>
44template <
int DisplacementDim>
56template <
int DisplacementDim>
60template <
int DisplacementDim>
92 :
value{1 + gamma_p * theta},
103template <
int DisplacementDim>
106 double const sqrtPhi,
double const alpha_p,
double const beta_p,
107 double const delta_p,
double const epsilon_p)
111 4 * boost::math::pow<2>(delta_p) * boost::math::pow<3>(s.
I_1)) /
113 3 * beta_p + 6 * epsilon_p * s.
I_1;
116template <
int DisplacementDim>
121 double const gamma_p,
double const m_p)
124 (s.
D + s.
J_2 * m_p * gamma_p * dtheta_dsigma / one_gt.
value)) /
128template <
int DisplacementDim>
133 double const I_1_squared = boost::math::pow<2>(s.
I_1);
136 return std::sqrt(s.
J_2 * std::pow(1 + mp.
gamma * s.
J_3 /
137 (s.
J_2 * std::sqrt(s.
J_2)),
139 mp.
alpha / 2. * I_1_squared +
140 boost::math::pow<2>(mp.
delta) *
141 boost::math::pow<2>(I_1_squared)) +
145template <
int DisplacementDim>
146std::pair<MathLib::KelvinVector::KelvinVectorType<DisplacementDim>,
151 constexpr int KelvinVectorSize =
161 return {KelvinVector::Zero(), KelvinMatrix::Zero()};
164 auto const& P_dev = Invariants::deviatoric_projection;
167 if (Invariants::determinant(s.
D) == 0)
169 OGS_FATAL(
"Determinant is zero. Matrix is non-invertable.");
173 KelvinVector
const sigma_D_inverse_D = P_dev * sigma_D_inverse;
175 KelvinVector
const dtheta_dsigma =
176 theta * sigma_D_inverse_D - 3. / 2. * theta / s.
J_2 * s.
D;
177 KelvinMatrix
const d2theta_dsigma2 =
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();
184 return {dtheta_dsigma, d2theta_dsigma2};
187template <
int DisplacementDim>
195 double const eps_p_V,
196 double const eps_p_V_dot,
197 double const eps_p_eff_dot,
202 static int const KelvinVectorSize =
208 auto const& identity2 = Invariants::identity2;
210 double const theta = s.
J_3 / (s.
J_2 * std::sqrt(s.
J_2));
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;
218 auto const [dtheta_dsigma, d2theta_dsigma2] =
222 double const sqrtPhi = std::sqrt(
224 boost::math::pow<2>(mp.
delta_p) * boost::math::pow<4>(s.
I_1));
226 s, one_gt, sqrtPhi, dtheta_dsigma, mp.
gamma_p, mp.
m_p);
227 KelvinVector
const lambda_flow_D = lambda * flow_D;
229 residual.template segment<KelvinVectorSize>(KelvinVectorSize).noalias() =
230 eps_p_D_dot - lambda_flow_D;
236 residual(2 * KelvinVectorSize, 0) = eps_p_V_dot - lambda * flow_V;
240 residual(2 * KelvinVectorSize + 1) =
242 std::sqrt(2. / 3. * lambda_flow_D.transpose() * lambda_flow_D);
245 residual(2 * KelvinVectorSize + 2) =
yieldFunction(mp, s, k) / mp.
G;
249template <
int DisplacementDim>
256 static int const KelvinVectorSize =
264 auto const& P_dev = Invariants::deviatoric_projection;
265 auto const& identity2 = Invariants::identity2;
267 double const theta = s.
J_3 / (s.
J_2 * std::sqrt(s.
J_2));
270 auto const [dtheta_dsigma, d2theta_dsigma2] =
274 double const sqrtPhi = std::sqrt(
276 boost::math::pow<2>(mp.
delta_p) * boost::math::pow<4>(s.
I_1));
278 s, one_gt, sqrtPhi, dtheta_dsigma, mp.
gamma_p, mp.
m_p);
279 KelvinVector
const lambda_flow_D = lambda * flow_D;
285 jacobian.template block<KelvinVectorSize, KelvinVectorSize>(0, 0)
286 .noalias() = KelvinMatrix::Identity();
290 .template block<KelvinVectorSize, KelvinVectorSize>(0, KelvinVectorSize)
291 .noalias() = 2 * KelvinMatrix::Identity();
294 jacobian.template block<KelvinVectorSize, 1>(0, 2 * KelvinVectorSize)
295 .noalias() = mp.
K / mp.
G * identity2;
303 KelvinVector
const M0 = s.
J_2 / one_gt.
value * dtheta_dsigma;
305 KelvinVector
const dPhi_dsigma =
306 one_gt.
pow_m_p * (s.
D + gm_p * M0) +
308 4 * boost::math::pow<2>(mp.
delta_p) * boost::math::pow<3>(s.
I_1)) *
312 KelvinMatrix
const M1 =
314 (s.
D * dPhi_dsigma.transpose() + gm_p * M0 * dPhi_dsigma.transpose());
316 KelvinMatrix
const M2 =
317 one_gt.
pow_m_p * (P_dev + s.
D * gm_p * M0.transpose());
320 KelvinMatrix
const M3 =
322 ((s.
D + (gm_p - mp.
gamma_p) * M0) * dtheta_dsigma.transpose() +
323 s.
J_2 * d2theta_dsigma2);
326 KelvinMatrix
const dflow_D_dsigma =
327 (-M1 / (4 * boost::math::pow<3>(sqrtPhi)) + (M2 + M3) / (2 * sqrtPhi)) *
330 .template block<KelvinVectorSize, KelvinVectorSize>(KelvinVectorSize, 0)
331 .noalias() = -lambda * dflow_D_dsigma;
335 .template block<KelvinVectorSize, KelvinVectorSize>(KelvinVectorSize,
337 .noalias() = KelvinMatrix::Identity() / dt;
343 .template block<KelvinVectorSize, 1>(KelvinVectorSize,
344 2 * KelvinVectorSize + 2)
345 .noalias() = -flow_D;
350 KelvinVector
const dflow_V_dsigma =
353 boost::math::pow<3>(s.
I_1)) /
354 (4 * boost::math::pow<3>(sqrtPhi)) * dPhi_dsigma +
356 12 * boost::math::pow<2>(mp.
delta_p * s.
I_1) * identity2) /
360 jacobian.template block<1, KelvinVectorSize>(2 * KelvinVectorSize, 0)
361 .noalias() = -lambda * dflow_V_dsigma.transpose();
367 jacobian(2 * KelvinVectorSize, 2 * KelvinVectorSize) = 1. / dt;
375 jacobian(2 * KelvinVectorSize, 2 * KelvinVectorSize + 2) = -flow_V;
379 double const eff_flow =
380 std::sqrt(2. / 3. * lambda_flow_D.transpose() * lambda_flow_D);
385 KelvinVector
const eff_flow23_lambda_flow_D =
386 -2 / 3. / eff_flow * lambda_flow_D;
389 .template block<1, KelvinVectorSize>(2 * KelvinVectorSize + 1, 0)
390 .noalias() = lambda * dflow_D_dsigma * eff_flow23_lambda_flow_D;
392 jacobian(2 * KelvinVectorSize + 1, 2 * KelvinVectorSize + 2) =
393 eff_flow23_lambda_flow_D.transpose() * flow_D;
399 jacobian(2 * KelvinVectorSize + 1, 2 * KelvinVectorSize + 1) = 1. / dt;
403 double const one_gt_pow_m = std::pow(one_gt.
value, mp.
m);
404 double const gm = mp.
gamma * mp.
m;
406 KelvinVector
const dF_dsigma =
408 (one_gt_pow_m * (s.
D + gm * M0) +
410 boost::math::pow<3>(s.
I_1)) *
416 .template block<1, KelvinVectorSize>(2 * KelvinVectorSize + 2, 0)
417 .noalias() = dF_dsigma.transpose() / mp.
G;
421 jacobian(2 * KelvinVectorSize + 2, 2 * KelvinVectorSize + 1) =
430template <
int DisplacementDim>
432 double const K,
double const G)
434 static int const KelvinVectorSize =
438 auto const& P_dev = Invariants::deviatoric_projection;
439 auto const& P_sph = Invariants::spherical_projection;
443 return -2. * I * P_dev - 3. * K / G * I * P_sph;
447 double const hardening_coefficient,
448 double const eps_p_eff)
450 return kappa * (1. + eps_p_eff * hardening_coefficient);
453template <
int DisplacementDim>
455 double const G,
double const K,
461 static int const KelvinVectorSize =
464 auto const& P_dev = Invariants::deviatoric_projection;
467 double const pressure_prev = Invariants::trace(sigma_prev) / (-3. * G);
469 double const e_prev = Invariants::trace(eps_prev);
471 double const pressure = pressure_prev - K / G * (eps_V - e_prev);
474 P_dev * sigma_prev / G;
477 sigma_D_prev + 2 * P_dev * (eps - eps_prev);
478 return sigma_D - pressure * Invariants::identity2;
483template <
typename Res
idualVector,
typename KelvinVector>
484std::tuple<KelvinVector, PlasticStrain<KelvinVector>,
double>
487 static auto const size = KelvinVector::SizeAtCompileTime;
488 return std::forward_as_tuple(
489 solution.template segment<size>(size * 0),
491 solution[size * 2], solution[size * 2 + 1]},
492 solution[size * 2 + 2]);
495template <
int DisplacementDim>
499 std::unique_ptr<DamagePropertiesParameters>&& damage_properties,
502 _mp(std::move(material_properties)),
508template <
int DisplacementDim>
516 material_state_variables)
const
519 &material_state_variables) !=
nullptr);
522 material_state_variables)
525 auto const& identity2 = Invariants::identity2;
526 return (eps - eps_p.D - eps_p.V / 3 * identity2).dot(sigma) / 2;
529template <
int DisplacementDim>
530std::optional<std::tuple<typename SolidEhlers<DisplacementDim>::KelvinVector,
532 DisplacementDim>::MaterialStateVariables>,
539 material_state_variables)
const
541 auto const& eps_m = std::get<MPL::SymmetricTensor<DisplacementDim>>(
543 auto const& eps_m_prev = std::get<MPL::SymmetricTensor<DisplacementDim>>(
545 auto const& sigma_prev = std::get<MPL::SymmetricTensor<DisplacementDim>>(
546 variable_array_prev.
stress);
549 &material_state_variables) !=
nullptr);
553 material_state_variables);
559 double const eps_V = Invariants::trace(eps_m);
561 auto const& P_dev = Invariants::deviatoric_projection;
569 mp.
G, mp.
K, sigma_prev, eps_m, eps_m_prev, eps_V);
579 if ((sigma.squaredNorm() == 0. || dt == 0. ||
583 state.
eps_p.eff)) < 0.))
586 mp.
K - 2. / 3 * mp.
G, mp.
G);
600 Eigen::Matrix<double, JacobianResidualSize, 1>;
614 auto const& eps_p_D =
615 solution.template segment<KelvinVectorSize>(
621 double const eps_p_V_dot = (eps_p_V - state.
eps_p_prev.V) / dt;
624 double const eps_p_eff_dot =
632 solution.template segment<KelvinVectorSize>(
645 auto const update_solution =
648 solution += increment;
650 mp.
G * solution.template segment<KelvinVectorSize>(0)};
654 linear_solver, update_jacobian, update_residual,
657 auto const success_iterations = newton_solver.solve(jacobian);
659 if (!success_iterations)
668 if (*success_iterations == 0)
670 linear_solver.compute(jacobian);
673 std::tie(sigma, state.
eps_p, std::ignore) =
682 Eigen::RowMajor>::Zero();
683 dresidual_deps.template block<KelvinVectorSize, KelvinVectorSize>(0, 0)
696 linear_solver.solve(-dresidual_deps)
697 .template block<KelvinVectorSize, KelvinVectorSize>(0, 0);
706 "Unimplemented tangent type behaviour for the tangent type "
714 return {std::make_tuple(
723template <
int DisplacementDim>
724std::vector<typename MechanicsBase<DisplacementDim>::InternalVariable>
727 return {{
"damage.kappa_d", 1,
730 std::vector<double>& cache) -> std::vector<double>
const&
734 auto const& ehlers_state =
742 state) -> std::span<double>
749 return {&ehlers_state.damage.kappa_d(), 1};
754 std::vector<double>& cache) -> std::vector<double>
const&
758 auto const& ehlers_state =
766 state) -> std::span<double>
773 return {&ehlers_state.damage.value(), 1};
775 {
"eps_p.D", KelvinVector::RowsAtCompileTime,
778 std::vector<double>& cache) -> std::vector<double>
const&
782 auto const& ehlers_state =
785 cache.resize(KelvinVector::RowsAtCompileTime);
787 cache, KelvinVector::RowsAtCompileTime) =
789 ehlers_state.eps_p.D);
794 state) -> std::span<double>
802 ehlers_state.eps_p.D.data(),
803 static_cast<std::size_t
>(KelvinVector::RowsAtCompileTime)};
808 std::vector<double>& cache) -> std::vector<double>
const&
812 auto const& ehlers_state =
816 cache.front() = ehlers_state.
eps_p.V;
820 state) -> std::span<double>
827 return {&ehlers_state.eps_p.V, 1};
832 std::vector<double>& cache) -> std::vector<double>
const&
836 auto const& ehlers_state =
840 cache.front() = ehlers_state.
eps_p.eff;
844 state) -> std::span<double>
851 return {&ehlers_state.eps_p.eff, 1};
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);
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.);
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);
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.;
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.;
895 result(5, 5) = v(0) * v(2) + v(5) * v(5) / 2.;
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);
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);
914 result(2, 2) = v(2) * v(2);
915 result(2, 3) = result(3, 2) = 0;
917 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.