OGS
LinearElasticIsotropicPhaseField.h
Go to the documentation of this file.
1
9
10#pragma once
11
12#include <Eigen/Eigenvalues>
13#include <boost/math/special_functions/pow.hpp>
14
16
17namespace MaterialLib
18{
19namespace Solids
20{
21namespace Phasefield
22{
29
30template <int DisplacementDim>
34
37inline double heaviside(double const v)
38{
39 return (v < 0) ? 0.0 : 1.0;
40}
41
44inline double macaulayTensile(double const v)
45{
46 return v * heaviside(v);
47}
48inline double macaulayCompressive(double v)
49{
50 return v * (1 - heaviside(v));
51}
52
53
55 int const i, int const j,
56 Eigen::Matrix<double, 3, 1> const& principal_strain);
57
59 int const i, int const j,
60 Eigen::Matrix<double, 3, 1> const& principal_strain);
61
62template <int DisplacementDim>
64 DisplacementDim> /* sigma_real */,
66 DisplacementDim> /* sigma_tensile */,
68 double /* strain_energy_tensile */, double /* elastic_energy */,
70 DisplacementDim> /* C_tensile */,
72 DisplacementDim> /* C_compressive
73 */
74 >
76 double const degradation, double const bulk_modulus, double const mu,
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}
138
139template <int DisplacementDim>
141 DisplacementDim> /* sigma_real */,
143 DisplacementDim> /* sigma_tensile */,
145 double /* strain_energy_tensile */, double /* elastic_energy */,
147 DisplacementDim> /* C_tensile */,
149 DisplacementDim> /* C_compressive
150 */
151 >
153 double const degradation,
154 double const lambda,
155 double const mu,
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}
255
256template <int DisplacementDim>
258 DisplacementDim> /* sigma_real */,
260 DisplacementDim> /* sigma_tensile */,
262 double /* strain_energy_tensile */, double /* elastic_energy */,
264 DisplacementDim> /* C_tensile */,
266 DisplacementDim> /* C_compressive
267 */
268 >
270 double const degradation,
271 double const bulk_modulus,
272 double const mu,
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}
308
309template <int DisplacementDim>
311 DisplacementDim> /* sigma_real */,
313 DisplacementDim> /* sigma_tensile */,
315 double /* strain_energy_tensile */, double /* elastic_energy */,
317 DisplacementDim> /* C_tensile */,
319 DisplacementDim> /* C_compressive
320 */
321 >
323 double const degradation,
324 double const bulk_modulus,
325 double const mu,
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}
376
377} // namespace Phasefield
378} // namespace Solids
379} // namespace MaterialLib
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)
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:
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)
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, 3, 3 > kelvinVectorToTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
KelvinVectorType< DisplacementDim > tensorToKelvin(Eigen::Matrix< double, 3, 3 > const &m)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor > KelvinMatrixType