70 double const degradation,
double const bulk_modulus,
double const mu,
73 static int const KelvinVectorSize =
81 auto const& P_dev = Invariants::deviatoric_projection;
82 KelvinVector
const epsd_curr = P_dev * eps;
85 double const eps_curr_trace = Invariants::trace(eps);
87 KelvinMatrix C_tensile = KelvinMatrix::Zero();
88 KelvinMatrix C_compressive = KelvinMatrix::Zero();
90 auto strain_energy_computation_vol = [&](
auto&& macaulay)
92 auto macaulay_squared = [&macaulay](
double x)
93 {
return boost::math::pow<2>(macaulay(x)); };
94 return bulk_modulus / 2 * macaulay_squared(eps_curr_trace);
97 auto stress_computation_vol = [&](
auto&& macaulay)
98 {
return bulk_modulus * macaulay(eps_curr_trace) * Invariants::identity2; };
100 auto hs = [&](
double const v) {
return heaviside(v); };
106 double const strain_energy_tensile = strain_energy_computation_vol(mt) +
107 mu * epsd_curr.transpose() * epsd_curr;
109 KelvinVector
const sigma_tensile =
110 stress_computation_vol(mt) + 2 * mu * epsd_curr;
112 KelvinVector
const sigma_compressive = stress_computation_vol(mc);
114 C_tensile.template topLeftCorner<3, 3>().setConstant(bulk_modulus *
116 C_tensile.noalias() += 2 * mu * P_dev * KelvinMatrix::Identity();
118 C_compressive.template topLeftCorner<3, 3>().setConstant(
119 bulk_modulus * (1 - hs(eps_curr_trace)));
121 double const elastic_energy =
122 bulk_modulus / 2 * eps_curr_trace * eps_curr_trace +
123 mu * epsd_curr.transpose() * epsd_curr;
124 KelvinVector
const sigma_real =
125 degradation * sigma_tensile + sigma_compressive;
127 KelvinMatrix
const D = degradation * C_tensile + C_compressive;
129 return std::make_tuple(sigma_real, sigma_tensile, D, strain_energy_tensile,
130 elastic_energy, C_tensile, C_compressive);
147 double const degradation,
152 static int const KelvinVectorSize =
160 KelvinMatrix C_tensile = KelvinMatrix::Zero();
161 KelvinMatrix C_compressive = KelvinMatrix::Zero();
166 Eigen::EigenSolver<
decltype(eps_tensor)> eigen_solver(eps_tensor);
167 Eigen::Matrix<double, 3, 1>
const principal_strain =
168 eigen_solver.eigenvalues().real();
169 double const sum_strain = principal_strain.sum();
171 std::array<KelvinVector, 3> M_kelvin;
173 for (
int i = 0; i < 3; i++)
176 eigen_solver.eigenvectors().real().col(i).normalized() *
177 eigen_solver.eigenvectors().real().col(i).normalized().transpose());
180 auto strain_energy_computation = [&](
auto&& macaulay)
182 auto macaulay_squared = [&macaulay](
double x)
183 {
return boost::math::pow<2>(macaulay(x)); };
184 return lambda / 2 * macaulay_squared(sum_strain) +
185 mu * principal_strain.unaryExpr(macaulay_squared).sum();
188 auto stress_computation = [&](
auto&& macaulay)
190 KelvinVector stress =
191 lambda * macaulay(sum_strain) * Invariants::identity2;
192 for (
int i = 0; i < 3; i++)
194 stress += 2 * mu * macaulay(principal_strain[i]) * M_kelvin[i];
199 auto hs = [&](
double const v) {
return heaviside(v); };
205 double const strain_energy_tensile = strain_energy_computation(mt);
207 double const strain_energy_compressive = strain_energy_computation(mc);
209 KelvinVector
const sigma_tensile = stress_computation(mt);
211 KelvinVector
const sigma_compressive = stress_computation(mc);
213 C_tensile.template topLeftCorner<3, 3>().setConstant(lambda *
216 for (
int i = 0; i < 3; i++)
218 C_tensile.noalias() += 2 * mu * hs(principal_strain[i]) * M_kelvin[i] *
219 M_kelvin[i].transpose();
220 for (
int j = 0; j < 3; j++)
222 C_tensile.noalias() +=
228 C_compressive.template topLeftCorner<3, 3>().setConstant(
229 lambda * (1 - hs(sum_strain)));
230 KelvinMatrix C_temp = KelvinMatrix::Zero();
231 for (
int i = 0; i < 3; i++)
233 C_compressive.noalias() += 2 * mu * (1 - hs(principal_strain[i])) *
234 M_kelvin[i] * M_kelvin[i].transpose();
235 C_temp.noalias() = M_kelvin[i] * M_kelvin[i].transpose();
236 for (
int j = 0; j < 3; j++)
238 C_compressive.noalias() +=
244 double const elastic_energy =
245 degradation * strain_energy_tensile + strain_energy_compressive;
247 KelvinVector sigma_real = degradation * sigma_tensile + sigma_compressive;
248 KelvinMatrix
const D = degradation * C_tensile + C_compressive;
250 return std::make_tuple(sigma_real, sigma_tensile, D, strain_energy_tensile,
251 elastic_energy, C_tensile, C_compressive);
268 double const degradation,
269 double const bulk_modulus,
273 static int const KelvinVectorSize =
281 auto const& P_dev = Invariants::deviatoric_projection;
282 KelvinVector
const epsd_curr = P_dev * eps;
285 double const eps_curr_trace = Invariants::trace(eps);
287 KelvinMatrix C_tensile = KelvinMatrix::Zero();
288 KelvinMatrix C_compressive = KelvinMatrix::Zero();
290 double const strain_energy_tensile =
291 bulk_modulus / 2 * eps_curr_trace * eps_curr_trace +
292 mu * epsd_curr.transpose() * epsd_curr;
293 double const elastic_energy = degradation * strain_energy_tensile;
294 KelvinVector
const sigma_tensile =
295 bulk_modulus * eps_curr_trace * Invariants::identity2 +
297 C_tensile.template topLeftCorner<3, 3>().setConstant(bulk_modulus);
298 C_tensile.noalias() += 2 * mu * P_dev * KelvinMatrix::Identity();
299 KelvinVector
const sigma_real = degradation * sigma_tensile;
301 KelvinMatrix
const D = degradation * C_tensile + C_compressive;
303 return std::make_tuple(sigma_real, sigma_tensile, D, strain_energy_tensile,
304 elastic_energy, C_tensile, C_compressive);
321 double const degradation,
322 double const bulk_modulus,
326 static int const KelvinVectorSize =
334 auto const& P_dev = Invariants::deviatoric_projection;
335 KelvinVector
const epsd_curr = P_dev * eps;
338 double const eps_curr_trace = Invariants::trace(eps);
340 KelvinMatrix C_tensile = KelvinMatrix::Zero();
341 KelvinMatrix C_compressive = KelvinMatrix::Zero();
343 KelvinVector
const sigma_tensile =
344 bulk_modulus * eps_curr_trace * Invariants::identity2 +
349 Eigen::SelfAdjointEigenSolver<
decltype(sigma_tensor)>
350 selfAdjoint_eigen_solver(sigma_tensor);
351 double const max_principle_stress =
352 selfAdjoint_eigen_solver
356 double const E0 = 9.0 * mu * bulk_modulus / (3.0 * bulk_modulus + mu);
357 auto hs = [&](
double const v) {
return heaviside(v); };
359 double const strain_energy_tensile = 0.5 * max_principle_stress *
360 max_principle_stress *
361 hs(max_principle_stress) / E0;
362 double const elastic_energy =
363 degradation * (bulk_modulus / 2 * eps_curr_trace * eps_curr_trace +
364 mu * epsd_curr.transpose() * epsd_curr);
366 C_tensile.template topLeftCorner<3, 3>().setConstant(bulk_modulus);
367 C_tensile.noalias() += 2 * mu * P_dev * KelvinMatrix::Identity();
368 KelvinVector
const sigma_real = degradation * sigma_tensile;
370 KelvinMatrix
const D = degradation * C_tensile + C_compressive;
371 return std::make_tuple(sigma_real, sigma_tensile, D, strain_energy_tensile,
372 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)