OGS
MaterialLib::Solids::Phasefield Namespace Reference

Classes

class  AT_DegradationDerivative
class  COHESIVE_DegradationDerivative
class  DegradationDerivative

Enumerations

enum class  PhaseFieldModel { AT1 , AT2 , COHESIVE }
enum class  SofteningCurve { Linear , Exponential }
enum class  EnergySplitModel {
  Isotropic , VolDev , Spectral , EffectiveStress ,
  OrthoVolDev , OrthoMasonry
}

Functions

double evaluateHTensSpectral (int const i, int const j, Eigen::Matrix< double, 3, 1 > const &principal_strain)
 H terms in the Spectral decomposition:
double evaluateHCompSpectral (int const i, int const j, Eigen::Matrix< double, 3, 1 > const &principal_strain)
template<>
MathLib::KelvinVector::KelvinMatrixType< 3 > aOdotB< 3 > (MathLib::KelvinVector::KelvinVectorType< 3 > const &A, MathLib::KelvinVector::KelvinVectorType< 3 > const &B)
template<>
MathLib::KelvinVector::KelvinMatrixType< 2 > aOdotB< 2 > (MathLib::KelvinVector::KelvinVectorType< 2 > const &A, MathLib::KelvinVector::KelvinVectorType< 2 > const &B)
template<int DisplacementDim>
MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > aOdotB (MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &A, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &B)
double heaviside (double const v)
double macaulayTensile (double const v)
double macaulayCompressive (double v)
template<int DisplacementDim>
std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > calculateVolDevDegradedStress (double const degradation, double const bulk_modulus, double const mu, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps)
template<int DisplacementDim>
std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > calculateSpectralDegradedStress (double const degradation, double const lambda, double const mu, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps)
template<int DisplacementDim>
std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > calculateIsotropicDegradedStress (double const degradation, double const bulk_modulus, double const mu, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps)
template<int DisplacementDim>
std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > calculateIsotropicDegradedStressWithRankineEnergy (double const degradation, double const bulk_modulus, double const mu, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps)
template<int DisplacementDim>
std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > calculateOrthoVolDevDegradedStress (double const degradation, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > const &C_ortho)
template<int DisplacementDim>
std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > calculateOrthoMasonryDegradedStress (double const degradation, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > const &C_ortho)
template<int DisplacementDim>
PhaseFieldModel convertStringToPhaseFieldModel (std::string const &phasefield_model)
template<int DisplacementDim>
SofteningCurve convertStringToSofteningCurve (std::optional< std::string > const &softening_curve)
template<int DisplacementDim>
EnergySplitModel convertStringToEnergySplitModel (std::string const &energy_split_model)
template<int DisplacementDim>
std::unique_ptr< DegradationDerivativecreatDegradationDerivative (PhaseFieldModel const &phasefield_model, double const lch, SofteningCurve const &softening_curve)
template<typename T_DNDX, typename T_N, typename T_W, typename T_D, typename T_LOCAL_JAC, typename T_LOCAL_RHS>
void calculateCrackLocalJacobianAndResidual (T_DNDX &dNdx, T_N &N, T_W &w, T_D &d, T_LOCAL_JAC &local_Jac, T_LOCAL_RHS local_rhs, double const gc, double const ls, PhaseFieldModel &phasefield_model)
template<typename T_VECTOR, typename T_MATRIX, int DisplacementDim>
void calculateStress (T_VECTOR &sigma, T_VECTOR &sigma_tensile, T_VECTOR &sigma_compressive, T_VECTOR &eps_tensile, T_MATRIX &D, T_MATRIX &C_tensile, T_MATRIX &C_compressive, double &strain_energy_tensile, double &elastic_energy, double const degradation, T_VECTOR const &eps, EnergySplitModel const &energy_split_model, double const t, ParameterLib::SpatialPosition const &x, MaterialLib::Solids::MechanicsBase< DisplacementDim > const &solid_material)

Enumeration Type Documentation

◆ EnergySplitModel

◆ PhaseFieldModel

◆ SofteningCurve

Function Documentation

◆ aOdotB()

template<int DisplacementDim>
MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > MaterialLib::Solids::Phasefield::aOdotB ( MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const & A,
MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const & B )

Decompose the stiffness into tensile and compressive part. Judging by the physical observations, compression perpendicular to a crack does not cause crack propagation. Thus, the phase-field parameter is only involved into the tensile part to degrade the elastic strain energy.

Referenced by calculateSpectralDegradedStress().

◆ aOdotB< 2 >()

Definition at line 33 of file LinearElasticIsotropicPhaseField.cpp.

119{
121
122 result(0, 0) = A(0) * B(0);
123 result(0, 1) = result(1, 0) = A(3) * B(3) / 2.;
124 result(0, 2) = result(2, 0) = 0;
125 result(0, 3) = result(3, 0) = (A(0) * B(3) + A(3) * B(0)) / 2;
126
127 result(1, 1) = A(1) * B(1);
128 result(1, 2) = result(2, 1) = 0;
129 result(1, 3) = result(3, 1) = (A(1) * B(3) + A(3) * B(1)) / 2.;
130
131 result(2, 2) = A(2) * B(2);
132 result(2, 3) = result(3, 2) = 0;
133
134 result(3, 3) = (A(0) * B(1) + A(1) * B(0) + A(3) * B(3)) / 2.;
135
136 return result;
137}
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor > KelvinMatrixType

References heaviside(), and macaulayCompressive().

◆ aOdotB< 3 >()

Double–minor symmetrized tensor product

aOdotB computes the symmetric 4th-order tensor product

\[(A \odot B)_{ijkl} = \tfrac{1}{4} (A_{ik}B_{jl} + A_{il}B_{jk} + A_{jk}B_{il} + A_{jl}B_{ik}), \]

which is symmetric in both index pairs \((i,j)\) and \((k,l)\) (“double-minor symmetry”).

Reference Mehrabadi, M.M., & Cowin, S.C. (1990). Eigentensors of linear anisotropic elastic materials. Q. J. Mech. Appl. Math., 43(1), 15–41.

Definition at line 33 of file LinearElasticIsotropicPhaseField.cpp.

74{
76
77 result(0, 0) = A(0) * B(0);
78 result(0, 1) = result(1, 0) = A(3) * B(3) / 2.;
79 result(0, 2) = result(2, 0) = A(5) * B(5) / 2.;
80 result(0, 3) = result(3, 0) = (A(0) * B(3) + A(3) * B(0)) / 2;
81 result(0, 4) = result(4, 0) =
82 (A(3) * B(5) + A(5) * B(3)) / (2. * std::sqrt(2.));
83 result(0, 5) = result(5, 0) = (A(0) * B(5) + A(5) * B(0)) / 2;
84
85 result(1, 1) = A(1) * B(1);
86 result(1, 2) = result(2, 1) = A(4) * B(4) / 2.;
87 result(1, 3) = result(3, 1) = (A(1) * B(3) + A(3) * B(1)) / 2.;
88 result(1, 4) = result(4, 1) = (A(1) * B(4) + A(4) * B(1)) / 2.;
89 result(1, 5) = result(5, 1) =
90 (A(3) * B(4) + A(4) * B(3)) / (2. * std::sqrt(2.));
91
92 result(2, 2) = A(2) * B(2);
93 result(2, 3) = result(3, 2) =
94 (A(4) * B(5) + A(5) * B(4)) / (2. * std::sqrt(2.));
95 result(2, 4) = result(4, 2) = (A(2) * B(4) + A(4) * B(2)) / 2.;
96 result(2, 5) = result(5, 2) = (A(2) * B(5) + A(5) * B(2)) / 2.;
97
98 result(3, 3) = (A(0) * B(1) + A(1) * B(0) + A(3) * B(3)) / 2.;
99 result(3, 4) = result(4, 3) =
100 (A(3) * B(4) + A(4) * B(3)) / 4. +
101 (A(1) * B(5) + A(5) * B(1)) / (2. * std::sqrt(2.));
102 result(3, 5) = result(5, 3) =
103 (A(3) * B(5) + A(5) * B(3)) / 4. +
104 (A(0) * B(4) + A(4) * B(0)) / (2. * std::sqrt(2.));
105
106 result(4, 4) = (A(1) * B(2) + A(2) * B(1) + A(4) * B(4)) / 2.;
107 result(4, 5) = result(5, 4) =
108 (A(4) * B(5) + A(5) * B(4)) / 4. +
109 (A(2) * B(3) + A(3) * B(2)) / (2. * std::sqrt(2.));
110
111 result(5, 5) = (A(0) * B(2) + A(2) * B(0) + A(5) * B(5)) / 2.;
112 return result;
113}

◆ calculateCrackLocalJacobianAndResidual()

template<typename T_DNDX, typename T_N, typename T_W, typename T_D, typename T_LOCAL_JAC, typename T_LOCAL_RHS>
void MaterialLib::Solids::Phasefield::calculateCrackLocalJacobianAndResidual ( T_DNDX & dNdx,
T_N & N,
T_W & w,
T_D & d,
T_LOCAL_JAC & local_Jac,
T_LOCAL_RHS local_rhs,
double const gc,
double const ls,
PhaseFieldModel & phasefield_model )

Definition at line 297 of file PhaseFieldBase.h.

302{
303 switch (phasefield_model)
304 {
306 {
307 auto const local_Jac_AT1 =
308 (gc * 0.75 * dNdx.transpose() * dNdx * ls * w).eval();
309 local_Jac.noalias() += local_Jac_AT1;
310
311 local_rhs.noalias() -=
312 gc * (-0.375 * N.transpose() / ls) * w + local_Jac_AT1 * d;
313 break;
314 }
316 {
317 auto const local_Jac_AT2 =
318 (gc * (N.transpose() * N / ls + dNdx.transpose() * dNdx * ls) *
319 w)
320 .eval();
321 local_Jac.noalias() += local_Jac_AT2;
322
323 local_rhs.noalias() -=
324 local_Jac_AT2 * d - gc * (N.transpose() / ls) * w;
325 break;
326 }
328 {
329 auto const local_Jac_COHESIVE =
330 (2.0 / std::numbers::pi * gc *
331 (-N.transpose() * N / ls + dNdx.transpose() * dNdx * ls) * w)
332 .eval();
333
334 local_Jac.noalias() += local_Jac_COHESIVE;
335
336 local_rhs.noalias() -= local_Jac_COHESIVE * d;
337 break;
338 }
339 default:
340 {
341 OGS_FATAL("Invalid phase field model specified.");
342 }
343 }
344};
#define OGS_FATAL(...)
Definition Error.h:26

References AT1, AT2, COHESIVE, and OGS_FATAL.

◆ calculateIsotropicDegradedStress()

template<int DisplacementDim>
std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > MaterialLib::Solids::Phasefield::calculateIsotropicDegradedStress ( double const degradation,
double const bulk_modulus,
double const mu,
MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const & eps )

Definition at line 269 of file LinearElasticIsotropicPhaseField.h.

274{
275 static int const KelvinVectorSize =
277 using KelvinVector =
279 using KelvinMatrix =
282 // calculation of deviatoric parts
283 auto const& P_dev = Invariants::deviatoric_projection;
284 KelvinVector const epsd_curr = P_dev * eps;
285
286 // Hydrostatic part for the stress and the tangent.
287 double const eps_curr_trace = Invariants::trace(eps);
288
289 KelvinMatrix C_tensile = KelvinMatrix::Zero();
290 KelvinMatrix C_compressive = KelvinMatrix::Zero();
291
292 double const strain_energy_tensile =
293 bulk_modulus / 2 * eps_curr_trace * eps_curr_trace +
294 mu * epsd_curr.transpose() * epsd_curr;
295 double const elastic_energy = degradation * strain_energy_tensile;
296 KelvinVector const sigma_tensile =
297 bulk_modulus * eps_curr_trace * Invariants::identity2 +
298 2 * mu * epsd_curr;
299 C_tensile.template topLeftCorner<3, 3>().setConstant(bulk_modulus);
300 C_tensile.noalias() += 2 * mu * P_dev * KelvinMatrix::Identity();
301 KelvinVector const sigma_real = degradation * sigma_tensile;
302
303 KelvinMatrix const D = degradation * C_tensile + C_compressive;
304
305 return std::make_tuple(sigma_real, sigma_tensile, D, strain_energy_tensile,
306 elastic_energy, C_tensile, C_compressive);
307}
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

References MathLib::KelvinVector::kelvin_vector_dimensions().

Referenced by calculateStress().

◆ calculateIsotropicDegradedStressWithRankineEnergy()

template<int DisplacementDim>
std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > MaterialLib::Solids::Phasefield::calculateIsotropicDegradedStressWithRankineEnergy ( double const degradation,
double const bulk_modulus,
double const mu,
MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const & eps )

Definition at line 322 of file LinearElasticIsotropicPhaseField.h.

327{
328 static int const KelvinVectorSize =
330 using KelvinVector =
332 using KelvinMatrix =
335 // calculation of deviatoric parts
336 auto const& P_dev = Invariants::deviatoric_projection;
337 KelvinVector const epsd_curr = P_dev * eps;
338
339 // Hydrostatic part for the stress and the tangent.
340 double const eps_curr_trace = Invariants::trace(eps);
341
342 KelvinMatrix C_tensile = KelvinMatrix::Zero();
343 KelvinMatrix C_compressive = KelvinMatrix::Zero();
344
345 KelvinVector const sigma_tensile =
346 bulk_modulus * eps_curr_trace * Invariants::identity2 +
347 2 * mu * epsd_curr;
348 // The maximum principal stress
349 auto sigma_tensor =
351 Eigen::SelfAdjointEigenSolver<decltype(sigma_tensor)>
352 selfAdjoint_eigen_solver(sigma_tensor);
353 double const max_principle_stress =
354 selfAdjoint_eigen_solver
355 .eigenvalues()[2]; // selfAdjoint_eigen_solver function gives an
356 // ascending sequence.
357 //
358 double const E0 = 9.0 * mu * bulk_modulus / (3.0 * bulk_modulus + mu);
359 auto hs = [&](double const v) { return heaviside(v); };
360 //
361 double const strain_energy_tensile = 0.5 * max_principle_stress *
362 max_principle_stress *
363 hs(max_principle_stress) / E0;
364 double const elastic_energy =
365 degradation * (bulk_modulus / 2 * eps_curr_trace * eps_curr_trace +
366 mu * epsd_curr.transpose() * epsd_curr);
367 //
368 C_tensile.template topLeftCorner<3, 3>().setConstant(bulk_modulus);
369 C_tensile.noalias() += 2 * mu * P_dev * KelvinMatrix::Identity();
370 KelvinVector const sigma_real = degradation * sigma_tensile;
371
372 KelvinMatrix const D = degradation * C_tensile + C_compressive;
373 return std::make_tuple(sigma_real, sigma_tensile, D, strain_energy_tensile,
374 elastic_energy, C_tensile, C_compressive);
375}
Eigen::Matrix< double, 3, 3 > kelvinVectorToTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)

References heaviside(), MathLib::KelvinVector::kelvin_vector_dimensions(), and MathLib::KelvinVector::kelvinVectorToTensor().

◆ calculateOrthoMasonryDegradedStress()

template<int DisplacementDim>
std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > MaterialLib::Solids::Phasefield::calculateOrthoMasonryDegradedStress ( double const degradation,
MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const & eps,
MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > const & C_ortho )

Definition at line 127 of file LinearElasticOrthotropicPhaseField.h.

131{
132 using KelvinVector =
134 using KelvinMatrix =
136
137 KelvinMatrix const C = C_ortho;
138
139 // calculation of square root of elasticity tensor
140 Eigen::SelfAdjointEigenSolver<KelvinMatrix> es(C);
141 KelvinMatrix const sqrt_C = es.operatorSqrt();
142 Eigen::SelfAdjointEigenSolver<KelvinMatrix> es_inverse(C.inverse());
143 KelvinMatrix const sqrt_C_inverse = es_inverse.operatorSqrt();
144
145 // strain in a transformed space
146 KelvinVector const epst = sqrt_C * eps;
147
148 // strain tensor in 3D
149 Eigen::Matrix3d const eps_3D =
151
152 Eigen::SelfAdjointEigenSolver<Eigen::Matrix3d> const es_eps_3D(eps_3D);
153 Eigen::Vector3d const eigen_values_eps_3D = es_eps_3D.eigenvalues().real();
154 Eigen::Matrix3d const eigen_vectors_eps_3D = es_eps_3D.eigenvectors();
155 Eigen::Matrix3d const E1_eigenp =
156 eigen_vectors_eps_3D.col(0) * eigen_vectors_eps_3D.col(0).transpose();
157
158 Eigen::Matrix3d const E3_eigenp =
159 eigen_vectors_eps_3D.col(2) * eigen_vectors_eps_3D.col(2).transpose();
160
161 KelvinVector const E1_vec =
163
164 KelvinVector const E3_vec =
166
167 KelvinMatrix const E1oE1 = E1_vec * E1_vec.transpose();
168 KelvinMatrix const E3oE3 = E3_vec * E3_vec.transpose();
169
170 KelvinMatrix I_p = KelvinMatrix::Zero();
171 KelvinMatrix I_n = KelvinMatrix::Zero();
172
173 KelvinMatrix const I_S = KelvinMatrix::Identity();
174 if (DisplacementDim == 2)
175 {
176 if (std::abs(eigen_values_eps_3D(2) - eigen_values_eps_3D(0)) >
177 std::numeric_limits<double>::epsilon())
178 {
179 I_p = (macaulayTensile(eigen_values_eps_3D(0)) -
180 macaulayTensile(eigen_values_eps_3D(2))) /
181 (eigen_values_eps_3D(0) - eigen_values_eps_3D(2)) *
182 (I_S - (E1oE1 + E3oE3)) +
183 (heaviside(eigen_values_eps_3D(0)) * E1oE1 +
184 heaviside(eigen_values_eps_3D(2)) * E3oE3);
185 I_n = (macaulayCompressive(eigen_values_eps_3D(0)) -
186 macaulayCompressive(eigen_values_eps_3D(2))) /
187 (eigen_values_eps_3D(0) - eigen_values_eps_3D(2)) *
188 (I_S - (E1oE1 + E3oE3)) +
189 (heaviside(-eigen_values_eps_3D(0)) * E1oE1 +
190 heaviside(-eigen_values_eps_3D(2)) * E3oE3);
191 }
192 else
193 {
194 I_p = heaviside(eigen_values_eps_3D(0)) * I_S;
195 I_n = heaviside(-eigen_values_eps_3D(0)) * I_S;
196 }
197 }
198
199 // strain tensile / compressive
200 KelvinVector const eps_tensile = sqrt_C_inverse * (I_p * sqrt_C) * eps;
201 KelvinVector const eps_compressive = sqrt_C_inverse * (I_n * sqrt_C) * eps;
202
203 // C tensile / compressive
204 KelvinMatrix const C_tensile = sqrt_C * I_p * sqrt_C;
205 KelvinMatrix const C_compressive = sqrt_C * I_n * sqrt_C;
206
207 // stress tensile / compressive
208 KelvinVector const sigma_tensile = C * eps_tensile;
209 KelvinVector const sigma_compressive = C * eps_compressive;
210
211 // decomposition of strain energy density
212 double const strain_energy_tensile =
213 0.5 * sigma_tensile.adjoint() * eps_tensile;
214 double const strain_energy_compressive =
215 0.5 * sigma_compressive.adjoint() * eps_compressive;
216
217 double const elastic_energy =
218 degradation * strain_energy_tensile + strain_energy_compressive;
219
220 KelvinVector const sigma_real =
221 degradation * sigma_tensile + sigma_compressive;
222
223 KelvinMatrix const D = degradation * C_tensile + C_compressive;
224
225 return std::make_tuple(eps_tensile, sigma_real, sigma_tensile, D,
226 strain_energy_tensile, elastic_energy, C_tensile,
227 C_compressive);
228}
KelvinVectorType< DisplacementDim > tensorToKelvin(Eigen::Matrix< double, 3, 3 > const &m)

References heaviside(), MathLib::KelvinVector::kelvinVectorToTensor(), macaulayCompressive(), macaulayTensile(), and MathLib::KelvinVector::tensorToKelvin().

Referenced by calculateStress().

◆ calculateOrthoVolDevDegradedStress()

template<int DisplacementDim>
std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > MaterialLib::Solids::Phasefield::calculateOrthoVolDevDegradedStress ( double const degradation,
MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const & eps,
MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > const & C_ortho )

Definition at line 38 of file LinearElasticOrthotropicPhaseField.h.

42{
43 static constexpr int KelvinVectorSize =
45 using KelvinVector =
47 using KelvinMatrix =
50 // calculation of deviatoric parts
51 auto const& P_dev = Invariants::deviatoric_projection;
52
53 KelvinMatrix const C = C_ortho;
54 // calculation of square root of elasticity tensor
55 Eigen::SelfAdjointEigenSolver<KelvinMatrix> es(C);
56 KelvinMatrix const sqrt_C = es.operatorSqrt();
57 Eigen::SelfAdjointEigenSolver<KelvinMatrix> es_inverse(C.inverse());
58 KelvinMatrix const sqrt_C_inverse = es_inverse.operatorSqrt();
59
60 // strain in a transformed space
61 KelvinVector const epst = sqrt_C * eps;
62 double const epst_curr_trace = Invariants::trace(epst);
63
64 // projection tensors in transformed space
65 KelvinMatrix teps_p = KelvinMatrix::Zero();
66 KelvinMatrix teps_n = KelvinMatrix::Zero();
67 if (epst_curr_trace >= 0) /* QQQ */
68 {
69 teps_p.template topLeftCorner<3, 3>().setConstant(1. / 3);
70 }
71 else
72 {
73 teps_n.template topLeftCorner<3, 3>().setConstant(1. / 3);
74 }
75
76 teps_p.noalias() += P_dev * KelvinMatrix::Identity();
77
78 // strain tensile / compressive
79 KelvinVector const eps_tensile = sqrt_C_inverse * (teps_p * sqrt_C) * eps;
80 KelvinVector const eps_compressive =
81 sqrt_C_inverse * (teps_n * sqrt_C) * eps;
82
83 // projection tensors in original space
84 KelvinMatrix const der_eps_p = sqrt_C_inverse * (teps_p * sqrt_C);
85 KelvinMatrix const der_eps_n = sqrt_C_inverse * (teps_n * sqrt_C);
86 // C tensile / compressive
87 KelvinMatrix const C_tensile = der_eps_p.transpose() * C * der_eps_p;
88 KelvinMatrix const C_compressive = der_eps_n.transpose() * C * der_eps_n;
89
90 // stress tensile / compressive
91 KelvinVector const sigma_tensile = C_tensile * eps;
92 KelvinVector const sigma_compressive = C_compressive * eps;
93
94 // decomposition of strain energy density
95 double const strain_energy_tensile =
96 0.5 * sigma_tensile.adjoint() * eps_tensile;
97 double const strain_energy_compressive =
98 0.5 * sigma_compressive.adjoint() * eps_compressive;
99
100 double const elastic_energy =
101 degradation * strain_energy_tensile + strain_energy_compressive;
102
103 KelvinVector const sigma_real =
104 degradation * sigma_tensile + sigma_compressive;
105 KelvinMatrix const D = degradation * C_tensile + C_compressive;
106
107 return std::make_tuple(eps_tensile, sigma_real, sigma_tensile,
108 sigma_compressive, D, strain_energy_tensile,
109 elastic_energy, C_tensile, C_compressive);
110}

References MathLib::KelvinVector::kelvin_vector_dimensions().

Referenced by calculateStress().

◆ calculateSpectralDegradedStress()

template<int DisplacementDim>
std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > MaterialLib::Solids::Phasefield::calculateSpectralDegradedStress ( double const degradation,
double const lambda,
double const mu,
MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const & eps )

Definition at line 152 of file LinearElasticIsotropicPhaseField.h.

157{
158 static int const KelvinVectorSize =
160 using KelvinVector =
162 using KelvinMatrix =
165
166 KelvinMatrix C_tensile = KelvinMatrix::Zero();
167 KelvinMatrix C_compressive = KelvinMatrix::Zero();
168
169 // non-const for Eigen solver.
170 auto eps_tensor = MathLib::KelvinVector::kelvinVectorToTensor(eps);
171
172 Eigen::EigenSolver<decltype(eps_tensor)> eigen_solver(eps_tensor);
173 Eigen::Matrix<double, 3, 1> const principal_strain =
174 eigen_solver.eigenvalues().real();
175 double const sum_strain = principal_strain.sum();
176
177 std::array<KelvinVector, 3> M_kelvin;
178
179 for (int i = 0; i < 3; i++)
180 {
182 eigen_solver.eigenvectors().real().col(i).normalized() *
183 eigen_solver.eigenvectors().real().col(i).normalized().transpose());
184 }
185
186 auto strain_energy_computation = [&](auto&& macaulay)
187 {
188 auto macaulay_squared = [&macaulay](double x)
189 { return boost::math::pow<2>(macaulay(x)); };
190 return lambda / 2 * macaulay_squared(sum_strain) +
191 mu * principal_strain.unaryExpr(macaulay_squared).sum();
192 };
193
194 auto stress_computation = [&](auto&& macaulay)
195 {
196 KelvinVector stress =
197 lambda * macaulay(sum_strain) * Invariants::identity2;
198 for (int i = 0; i < 3; i++)
199 stress += 2 * mu * macaulay(principal_strain[i]) * M_kelvin[i];
200 return stress;
201 };
202
203 auto hs = [&](double const v) { return heaviside(v); };
204
205 auto mt = [&](double const v) { return macaulayTensile(v); };
206
207 auto mc = [&](double const v) { return macaulayCompressive(v); };
208
209 double const strain_energy_tensile = strain_energy_computation(mt);
210
211 double const strain_energy_compressive = strain_energy_computation(mc);
212
213 KelvinVector const sigma_tensile = stress_computation(mt);
214
215 KelvinVector const sigma_compressive = stress_computation(mc);
216
217 C_tensile.template topLeftCorner<3, 3>().setConstant(lambda *
218 hs(sum_strain));
219
220 for (int i = 0; i < 3; i++)
221 {
222 C_tensile.noalias() += 2 * mu * hs(principal_strain[i]) * M_kelvin[i] *
223 M_kelvin[i].transpose();
224 for (int j = 0; j < 3; j++)
225 {
226 C_tensile.noalias() +=
227 mu * evaluateHTensSpectral(i, j, principal_strain) *
228 aOdotB<DisplacementDim>(M_kelvin[i], M_kelvin[j]);
229 }
230 }
231
232 C_compressive.template topLeftCorner<3, 3>().setConstant(
233 lambda * (1 - hs(sum_strain)));
234 KelvinMatrix C_temp = KelvinMatrix::Zero();
235 for (int i = 0; i < 3; i++)
236 {
237 C_compressive.noalias() += 2 * mu * (1 - hs(principal_strain[i])) *
238 M_kelvin[i] * M_kelvin[i].transpose();
239 C_temp.noalias() = M_kelvin[i] * M_kelvin[i].transpose();
240 for (int j = 0; j < 3; j++)
241 C_compressive.noalias() +=
242 mu * evaluateHCompSpectral(i, j, principal_strain) *
243 aOdotB<DisplacementDim>(M_kelvin[i], M_kelvin[j]);
244 }
245
246 double const elastic_energy =
247 degradation * strain_energy_tensile + strain_energy_compressive;
248
249 KelvinVector sigma_real = degradation * sigma_tensile + sigma_compressive;
250 KelvinMatrix const D = degradation * C_tensile + C_compressive;
251
252 return std::make_tuple(sigma_real, sigma_tensile, D, strain_energy_tensile,
253 elastic_energy, C_tensile, C_compressive);
254}
double evaluateHCompSpectral(int const i, int const j, Eigen::Matrix< double, 3, 1 > const &principal_strain)
MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > aOdotB(MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &A, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &B)
double evaluateHTensSpectral(int const i, int const j, Eigen::Matrix< double, 3, 1 > const &principal_strain)
H terms in the Spectral decomposition:
KV::KelvinMatrixType< DisplacementDim > KelvinMatrix
KV::KelvinVectorType< DisplacementDim > KelvinVector

References aOdotB(), evaluateHCompSpectral(), evaluateHTensSpectral(), heaviside(), MathLib::KelvinVector::kelvin_vector_dimensions(), MathLib::KelvinVector::kelvinVectorToTensor(), macaulayCompressive(), macaulayTensile(), and MathLib::KelvinVector::tensorToKelvin().

Referenced by calculateStress().

◆ calculateStress()

template<typename T_VECTOR, typename T_MATRIX, int DisplacementDim>
void MaterialLib::Solids::Phasefield::calculateStress ( T_VECTOR & sigma,
T_VECTOR & sigma_tensile,
T_VECTOR & sigma_compressive,
T_VECTOR & eps_tensile,
T_MATRIX & D,
T_MATRIX & C_tensile,
T_MATRIX & C_compressive,
double & strain_energy_tensile,
double & elastic_energy,
double const degradation,
T_VECTOR const & eps,
EnergySplitModel const & energy_split_model,
double const t,
ParameterLib::SpatialPosition const & x,
MaterialLib::Solids::MechanicsBase< DisplacementDim > const & solid_material )

Definition at line 347 of file PhaseFieldBase.h.

355{
356 auto linear_elastic_mp =
358 DisplacementDim> const&>(solid_material)
359 .getMaterialProperties();
360
361 auto const lambda = linear_elastic_mp.lambda(t, x);
362 auto const bulk_modulus = linear_elastic_mp.bulk_modulus(t, x);
363 auto const mu = linear_elastic_mp.mu(t, x);
364 switch (energy_split_model)
365 {
367 {
368 std::tie(sigma, sigma_tensile, D, strain_energy_tensile,
369 elastic_energy, C_tensile, C_compressive) =
372 degradation, bulk_modulus, mu, eps);
373 break;
374 }
376 {
377 std::tie(sigma, sigma_tensile, D, strain_energy_tensile,
378 elastic_energy, C_tensile, C_compressive) =
380 DisplacementDim>(degradation, bulk_modulus, mu, eps);
381 break;
382 }
384 {
385 std::tie(sigma, sigma_tensile, D, strain_energy_tensile,
386 elastic_energy, C_tensile, C_compressive) =
389 degradation, lambda, mu, eps);
390 break;
391 }
393 {
394 std::tie(sigma, sigma_tensile, D, strain_energy_tensile,
395 elastic_energy, C_tensile, C_compressive) =
396 MaterialLib::Solids::Phasefield::
397 calculateIsotropicDegradedStressWithRankineEnergy<
398 DisplacementDim>(degradation, bulk_modulus, mu, eps);
399 break;
400 }
402 {
403 double temperature = 0.;
405 static_cast<MaterialLib::Solids::LinearElasticOrthotropic<
406 DisplacementDim> const&>(solid_material)
407 .getElasticTensor(t, x, temperature);
408
409 std::tie(eps_tensile, sigma, sigma_tensile, sigma_compressive, D,
410 strain_energy_tensile, elastic_energy, C_tensile,
411 C_compressive) = MaterialLib::Solids::Phasefield::
413 degradation, eps, C_ortho);
414 break;
415 }
417 {
418 double temperature = 0.;
420 static_cast<MaterialLib::Solids::LinearElasticOrthotropic<
421 DisplacementDim> const&>(solid_material)
422 .getElasticTensor(t, x, temperature);
423
424 std::tie(eps_tensile, sigma, sigma_tensile, D,
425 strain_energy_tensile, elastic_energy, C_tensile,
426 C_compressive) = MaterialLib::Solids::Phasefield::
428 degradation, eps, C_ortho);
429 break;
430 }
431 default:
432 OGS_FATAL("Invalid energy split model specified.");
433 }
434};
std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > calculateOrthoVolDevDegradedStress(double const degradation, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > const &C_ortho)
std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > calculateIsotropicDegradedStress(double const degradation, double const bulk_modulus, double const mu, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps)
std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > calculateSpectralDegradedStress(double const degradation, double const lambda, double const mu, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps)
std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > calculateVolDevDegradedStress(double const degradation, double const bulk_modulus, double const mu, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps)
std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > calculateOrthoMasonryDegradedStress(double const degradation, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > const &C_ortho)

References calculateIsotropicDegradedStress(), calculateOrthoMasonryDegradedStress(), calculateOrthoVolDevDegradedStress(), calculateSpectralDegradedStress(), calculateVolDevDegradedStress(), EffectiveStress, Isotropic, OGS_FATAL, OrthoMasonry, OrthoVolDev, Spectral, and VolDev.

Referenced by ProcessLib::HMPhaseField::IntegrationPointData< BMatricesType, ShapeMatricesType, DisplacementDim >::updateConstitutiveRelation(), and ProcessLib::PhaseField::IntegrationPointData< BMatricesType, ShapeMatricesType, DisplacementDim >::updateConstitutiveRelation().

◆ calculateVolDevDegradedStress()

template<int DisplacementDim>
std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > MaterialLib::Solids::Phasefield::calculateVolDevDegradedStress ( double const degradation,
double const bulk_modulus,
double const mu,
MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const & eps )

Definition at line 75 of file LinearElasticIsotropicPhaseField.h.

78{
79 static int const KelvinVectorSize =
81 using KelvinVector =
83 using KelvinMatrix =
86 // calculation of deviatoric parts
87 auto const& P_dev = Invariants::deviatoric_projection;
88 KelvinVector const epsd_curr = P_dev * eps;
89
90 // Hydrostatic part for the stress and the tangent.
91 double const eps_curr_trace = Invariants::trace(eps);
92
93 KelvinMatrix C_tensile = KelvinMatrix::Zero();
94 KelvinMatrix C_compressive = KelvinMatrix::Zero();
95
96 auto strain_energy_computation_vol = [&](auto&& macaulay)
97 {
98 auto macaulay_squared = [&macaulay](double x)
99 { return boost::math::pow<2>(macaulay(x)); };
100 return bulk_modulus / 2 * macaulay_squared(eps_curr_trace);
101 };
102
103 auto stress_computation_vol = [&](auto&& macaulay)
104 { return bulk_modulus * macaulay(eps_curr_trace) * Invariants::identity2; };
105
106 auto hs = [&](double const v) { return heaviside(v); };
107
108 auto mt = [&](double const v) { return macaulayTensile(v); };
109
110 auto mc = [&](double const v) { return macaulayCompressive(v); };
111
112 double const strain_energy_tensile = strain_energy_computation_vol(mt) +
113 mu * epsd_curr.transpose() * epsd_curr;
114
115 KelvinVector const sigma_tensile =
116 stress_computation_vol(mt) + 2 * mu * epsd_curr;
117
118 KelvinVector const sigma_compressive = stress_computation_vol(mc);
119
120 C_tensile.template topLeftCorner<3, 3>().setConstant(bulk_modulus *
121 hs(eps_curr_trace));
122 C_tensile.noalias() += 2 * mu * P_dev * KelvinMatrix::Identity();
123
124 C_compressive.template topLeftCorner<3, 3>().setConstant(
125 bulk_modulus * (1 - hs(eps_curr_trace)));
126
127 double const elastic_energy =
128 bulk_modulus / 2 * eps_curr_trace * eps_curr_trace +
129 mu * epsd_curr.transpose() * epsd_curr;
130 KelvinVector const sigma_real =
131 degradation * sigma_tensile + sigma_compressive;
132
133 KelvinMatrix const D = degradation * C_tensile + C_compressive;
134
135 return std::make_tuple(sigma_real, sigma_tensile, D, strain_energy_tensile,
136 elastic_energy, C_tensile, C_compressive);
137}

References heaviside(), MathLib::KelvinVector::kelvin_vector_dimensions(), macaulayCompressive(), and macaulayTensile().

Referenced by calculateStress().

◆ convertStringToEnergySplitModel()

template<int DisplacementDim>
EnergySplitModel MaterialLib::Solids::Phasefield::convertStringToEnergySplitModel ( std::string const & energy_split_model)

Definition at line 96 of file PhaseFieldBase.h.

98{
99 if (energy_split_model == "Isotropic")
100 {
102 }
103 if (energy_split_model == "VolumetricDeviatoric")
104 {
106 }
107 if (energy_split_model == "Spectral")
108 {
110 }
111 if (energy_split_model == "EffectiveStress")
112 {
114 }
115 if (energy_split_model == "OrthoVolDev")
116 {
118 }
119 if (energy_split_model == "OrthoMasonry")
120 {
122 }
123 OGS_FATAL(
124 "energy_split_model must be 'Isotropic', 'VolumetricDeviatoric', "
125 "'EffectiveStress', 'OrthoVolDev' or 'OrthoMasonry' but '{}' was given",
126 energy_split_model.c_str());
127};

References EffectiveStress, Isotropic, OGS_FATAL, OrthoMasonry, OrthoVolDev, Spectral, and VolDev.

Referenced by ProcessLib::HMPhaseField::createHMPhaseFieldProcess(), and ProcessLib::PhaseField::createPhaseFieldProcess().

◆ convertStringToPhaseFieldModel()

template<int DisplacementDim>
PhaseFieldModel MaterialLib::Solids::Phasefield::convertStringToPhaseFieldModel ( std::string const & phasefield_model)

Definition at line 36 of file PhaseFieldBase.h.

38{
39 if (phasefield_model == "AT1")
40 {
42 }
43 if (phasefield_model == "AT2")
44 {
46 }
47 if (phasefield_model == "COHESIVE")
48 {
50 }
52 "phasefield_model must be 'AT1', 'AT2' or 'COHESIVE' but '{}' "
53 "was given",
54 phasefield_model.c_str());
55};

References AT1, AT2, COHESIVE, and OGS_FATAL.

Referenced by ProcessLib::HMPhaseField::createHMPhaseFieldProcess(), and ProcessLib::PhaseField::createPhaseFieldProcess().

◆ convertStringToSofteningCurve()

template<int DisplacementDim>
SofteningCurve MaterialLib::Solids::Phasefield::convertStringToSofteningCurve ( std::optional< std::string > const & softening_curve)

Definition at line 64 of file PhaseFieldBase.h.

66{
67 if (softening_curve)
68 {
69 if (*softening_curve == "Linear")
70 {
72 }
73 if (*softening_curve == "Exponential")
74 {
76 }
78 "softening_curve must be 'Linear' or 'Exponential' but '{}' "
79 "was given",
80 softening_curve->c_str());
81 }
82 return SofteningCurve::Linear; // default
83};

References Exponential, Linear, and OGS_FATAL.

Referenced by ProcessLib::HMPhaseField::createHMPhaseFieldProcess(), and ProcessLib::PhaseField::createPhaseFieldProcess().

◆ creatDegradationDerivative()

template<int DisplacementDim>
std::unique_ptr< DegradationDerivative > MaterialLib::Solids::Phasefield::creatDegradationDerivative ( PhaseFieldModel const & phasefield_model,
double const lch,
SofteningCurve const & softening_curve )

Definition at line 281 of file PhaseFieldBase.h.

284{
285 switch (phasefield_model)
286 {
288 return std::make_unique<COHESIVE_DegradationDerivative>(
289 lch, softening_curve);
290 default:
291 return std::make_unique<AT_DegradationDerivative>();
292 }
293};

References COHESIVE.

◆ evaluateHCompSpectral()

double MaterialLib::Solids::Phasefield::evaluateHCompSpectral ( int const i,
int const j,
Eigen::Matrix< double, 3, 1 > const & principal_strain )

Definition at line 33 of file LinearElasticIsotropicPhaseField.cpp.

36{
37 if (i == j)
38 {
39 return 0.0;
40 }
41 double num_zero = 1.e-14;
42
43 if (std::fabs(principal_strain[i] - principal_strain[j]) < num_zero)
44 {
45 return 2 * (1 - heaviside(principal_strain[i]));
46 }
47 return 2 *
48 (macaulayCompressive(principal_strain[i]) -
49 macaulayCompressive(principal_strain[j])) /
50 (principal_strain[i] - principal_strain[j]);
51}

Referenced by calculateSpectralDegradedStress().

◆ evaluateHTensSpectral()

double MaterialLib::Solids::Phasefield::evaluateHTensSpectral ( int const i,
int const j,
Eigen::Matrix< double, 3, 1 > const & principal_strain )

H terms in the Spectral decomposition:

Definition at line 14 of file LinearElasticIsotropicPhaseField.cpp.

17{
18 if (i == j)
19 {
20 return 0.0;
21 }
22 if (std::fabs(principal_strain[i] - principal_strain[j]) <
23 std::numeric_limits<double>::epsilon())
24 {
25 return 2 * heaviside(principal_strain[i]);
26 }
27 return 2 *
28 (macaulayTensile(principal_strain[i]) -
29 macaulayTensile(principal_strain[j])) /
30 (principal_strain[i] - principal_strain[j]);
31}

References heaviside(), and macaulayTensile().

Referenced by calculateSpectralDegradedStress().

◆ heaviside()

double MaterialLib::Solids::Phasefield::heaviside ( double const v)
inline

heaviside function returns 1.0 if the argument is positive and 0.0 if negative

Definition at line 37 of file LinearElasticIsotropicPhaseField.h.

38{
39 return (v < 0) ? 0.0 : 1.0;
40}

Referenced by aOdotB< 2 >(), calculateIsotropicDegradedStressWithRankineEnergy(), calculateOrthoMasonryDegradedStress(), calculateSpectralDegradedStress(), calculateVolDevDegradedStress(), evaluateHTensSpectral(), macaulayCompressive(), and macaulayTensile().

◆ macaulayCompressive()

double MaterialLib::Solids::Phasefield::macaulayCompressive ( double v)
inline

◆ macaulayTensile()

double MaterialLib::Solids::Phasefield::macaulayTensile ( double const v)
inline

Macaulay brackets: positive strain is tensile and negative strain for compressive

Definition at line 44 of file LinearElasticIsotropicPhaseField.h.

45{
46 return v * heaviside(v);
47}

References heaviside().

Referenced by calculateOrthoMasonryDegradedStress(), calculateSpectralDegradedStress(), calculateVolDevDegradedStress(), and evaluateHTensSpectral().