76 double const degradation,
double const bulk_modulus,
double const mu,
79 static int const KelvinVectorSize =
87 auto const& P_dev = Invariants::deviatoric_projection;
88 KelvinVector
const epsd_curr = P_dev * eps;
91 double const eps_curr_trace = Invariants::trace(eps);
93 KelvinMatrix C_tensile = KelvinMatrix::Zero();
94 KelvinMatrix C_compressive = KelvinMatrix::Zero();
96 auto strain_energy_computation_vol = [&](
auto&& macaulay)
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);
103 auto stress_computation_vol = [&](
auto&& macaulay)
104 {
return bulk_modulus * macaulay(eps_curr_trace) * Invariants::identity2; };
106 auto hs = [&](
double const v) {
return heaviside(v); };
112 double const strain_energy_tensile = strain_energy_computation_vol(mt) +
113 mu * epsd_curr.transpose() * epsd_curr;
115 KelvinVector
const sigma_tensile =
116 stress_computation_vol(mt) + 2 * mu * epsd_curr;
118 KelvinVector
const sigma_compressive = stress_computation_vol(mc);
120 C_tensile.template topLeftCorner<3, 3>().setConstant(bulk_modulus *
122 C_tensile.noalias() += 2 * mu * P_dev * KelvinMatrix::Identity();
124 C_compressive.template topLeftCorner<3, 3>().setConstant(
125 bulk_modulus * (1 - hs(eps_curr_trace)));
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;
133 KelvinMatrix
const D = degradation * C_tensile + C_compressive;
135 return std::make_tuple(sigma_real, sigma_tensile, D, strain_energy_tensile,
136 elastic_energy, C_tensile, C_compressive);
153 double const degradation,
158 static int const KelvinVectorSize =
166 KelvinMatrix C_tensile = KelvinMatrix::Zero();
167 KelvinMatrix C_compressive = KelvinMatrix::Zero();
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();
177 std::array<KelvinVector, 3> M_kelvin;
179 for (
int i = 0; i < 3; i++)
182 eigen_solver.eigenvectors().real().col(i).normalized() *
183 eigen_solver.eigenvectors().real().col(i).normalized().transpose());
186 auto strain_energy_computation = [&](
auto&& macaulay)
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();
194 auto stress_computation = [&](
auto&& macaulay)
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];
203 auto hs = [&](
double const v) {
return heaviside(v); };
209 double const strain_energy_tensile = strain_energy_computation(mt);
211 double const strain_energy_compressive = strain_energy_computation(mc);
213 KelvinVector
const sigma_tensile = stress_computation(mt);
215 KelvinVector
const sigma_compressive = stress_computation(mc);
217 C_tensile.template topLeftCorner<3, 3>().setConstant(lambda *
220 for (
int i = 0; i < 3; i++)
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++)
226 C_tensile.noalias() +=
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++)
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() +=
246 double const elastic_energy =
247 degradation * strain_energy_tensile + strain_energy_compressive;
249 KelvinVector sigma_real = degradation * sigma_tensile + sigma_compressive;
250 KelvinMatrix
const D = degradation * C_tensile + C_compressive;
252 return std::make_tuple(sigma_real, sigma_tensile, D, strain_energy_tensile,
253 elastic_energy, C_tensile, C_compressive);
270 double const degradation,
271 double const bulk_modulus,
275 static int const KelvinVectorSize =
283 auto const& P_dev = Invariants::deviatoric_projection;
284 KelvinVector
const epsd_curr = P_dev * eps;
287 double const eps_curr_trace = Invariants::trace(eps);
289 KelvinMatrix C_tensile = KelvinMatrix::Zero();
290 KelvinMatrix C_compressive = KelvinMatrix::Zero();
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 +
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;
303 KelvinMatrix
const D = degradation * C_tensile + C_compressive;
305 return std::make_tuple(sigma_real, sigma_tensile, D, strain_energy_tensile,
306 elastic_energy, C_tensile, C_compressive);
323 double const degradation,
324 double const bulk_modulus,
328 static int const KelvinVectorSize =
336 auto const& P_dev = Invariants::deviatoric_projection;
337 KelvinVector
const epsd_curr = P_dev * eps;
340 double const eps_curr_trace = Invariants::trace(eps);
342 KelvinMatrix C_tensile = KelvinMatrix::Zero();
343 KelvinMatrix C_compressive = KelvinMatrix::Zero();
345 KelvinVector
const sigma_tensile =
346 bulk_modulus * eps_curr_trace * Invariants::identity2 +
351 Eigen::SelfAdjointEigenSolver<
decltype(sigma_tensor)>
352 selfAdjoint_eigen_solver(sigma_tensor);
353 double const max_principle_stress =
354 selfAdjoint_eigen_solver
358 double const E0 = 9.0 * mu * bulk_modulus / (3.0 * bulk_modulus + mu);
359 auto hs = [&](
double const v) {
return heaviside(v); };
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);
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;
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);
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::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)