OGS
LinearElasticIsotropicPhaseField.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
6#include <Eigen/Eigenvalues>
7#include <boost/math/special_functions/pow.hpp>
8
10
11namespace MaterialLib
12{
13namespace Solids
14{
15namespace Phasefield
16{
23
24template <int DisplacementDim>
28
31inline double heaviside(double const v)
32{
33 return (v < 0) ? 0.0 : 1.0;
34}
35
38inline double macaulayTensile(double const v)
39{
40 return v * heaviside(v);
41}
42inline double macaulayCompressive(double v)
43{
44 return v * (1 - heaviside(v));
45}
46
47
49 int const i, int const j,
50 Eigen::Matrix<double, 3, 1> const& principal_strain);
51
53 int const i, int const j,
54 Eigen::Matrix<double, 3, 1> const& principal_strain);
55
56template <int DisplacementDim>
58 DisplacementDim> /* sigma_real */,
60 DisplacementDim> /* sigma_tensile */,
62 double /* strain_energy_tensile */, double /* elastic_energy */,
64 DisplacementDim> /* C_tensile */,
66 DisplacementDim> /* C_compressive
67 */
68 >
70 double const degradation, double const bulk_modulus, double const mu,
72{
73 static int const KelvinVectorSize =
75 using KelvinVector =
77 using KelvinMatrix =
80 // calculation of deviatoric parts
81 auto const& P_dev = Invariants::deviatoric_projection;
82 KelvinVector const epsd_curr = P_dev * eps;
83
84 // Hydrostatic part for the stress and the tangent.
85 double const eps_curr_trace = Invariants::trace(eps);
86
87 KelvinMatrix C_tensile = KelvinMatrix::Zero();
88 KelvinMatrix C_compressive = KelvinMatrix::Zero();
89
90 auto strain_energy_computation_vol = [&](auto&& macaulay)
91 {
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);
95 };
96
97 auto stress_computation_vol = [&](auto&& macaulay)
98 { return bulk_modulus * macaulay(eps_curr_trace) * Invariants::identity2; };
99
100 auto hs = [&](double const v) { return heaviside(v); };
101
102 auto mt = [&](double const v) { return macaulayTensile(v); };
103
104 auto mc = [&](double const v) { return macaulayCompressive(v); };
105
106 double const strain_energy_tensile = strain_energy_computation_vol(mt) +
107 mu * epsd_curr.transpose() * epsd_curr;
108
109 KelvinVector const sigma_tensile =
110 stress_computation_vol(mt) + 2 * mu * epsd_curr;
111
112 KelvinVector const sigma_compressive = stress_computation_vol(mc);
113
114 C_tensile.template topLeftCorner<3, 3>().setConstant(bulk_modulus *
115 hs(eps_curr_trace));
116 C_tensile.noalias() += 2 * mu * P_dev * KelvinMatrix::Identity();
117
118 C_compressive.template topLeftCorner<3, 3>().setConstant(
119 bulk_modulus * (1 - hs(eps_curr_trace)));
120
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;
126
127 KelvinMatrix const D = degradation * C_tensile + C_compressive;
128
129 return std::make_tuple(sigma_real, sigma_tensile, D, strain_energy_tensile,
130 elastic_energy, C_tensile, C_compressive);
131}
132
133template <int DisplacementDim>
135 DisplacementDim> /* sigma_real */,
137 DisplacementDim> /* sigma_tensile */,
139 double /* strain_energy_tensile */, double /* elastic_energy */,
141 DisplacementDim> /* C_tensile */,
143 DisplacementDim> /* C_compressive
144 */
145 >
147 double const degradation,
148 double const lambda,
149 double const mu,
151{
152 static int const KelvinVectorSize =
154 using KelvinVector =
156 using KelvinMatrix =
159
160 KelvinMatrix C_tensile = KelvinMatrix::Zero();
161 KelvinMatrix C_compressive = KelvinMatrix::Zero();
162
163 // non-const for Eigen solver.
164 auto eps_tensor = MathLib::KelvinVector::kelvinVectorToTensor(eps);
165
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();
170
171 std::array<KelvinVector, 3> M_kelvin;
172
173 for (int i = 0; i < 3; i++)
174 {
176 eigen_solver.eigenvectors().real().col(i).normalized() *
177 eigen_solver.eigenvectors().real().col(i).normalized().transpose());
178 }
179
180 auto strain_energy_computation = [&](auto&& macaulay)
181 {
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();
186 };
187
188 auto stress_computation = [&](auto&& macaulay)
189 {
190 KelvinVector stress =
191 lambda * macaulay(sum_strain) * Invariants::identity2;
192 for (int i = 0; i < 3; i++)
193 {
194 stress += 2 * mu * macaulay(principal_strain[i]) * M_kelvin[i];
195 }
196 return stress;
197 };
198
199 auto hs = [&](double const v) { return heaviside(v); };
200
201 auto mt = [&](double const v) { return macaulayTensile(v); };
202
203 auto mc = [&](double const v) { return macaulayCompressive(v); };
204
205 double const strain_energy_tensile = strain_energy_computation(mt);
206
207 double const strain_energy_compressive = strain_energy_computation(mc);
208
209 KelvinVector const sigma_tensile = stress_computation(mt);
210
211 KelvinVector const sigma_compressive = stress_computation(mc);
212
213 C_tensile.template topLeftCorner<3, 3>().setConstant(lambda *
214 hs(sum_strain));
215
216 for (int i = 0; i < 3; i++)
217 {
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++)
221 {
222 C_tensile.noalias() +=
223 mu * evaluateHTensSpectral(i, j, principal_strain) *
224 aOdotB<DisplacementDim>(M_kelvin[i], M_kelvin[j]);
225 }
226 }
227
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++)
232 {
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++)
237 {
238 C_compressive.noalias() +=
239 mu * evaluateHCompSpectral(i, j, principal_strain) *
240 aOdotB<DisplacementDim>(M_kelvin[i], M_kelvin[j]);
241 }
242 }
243
244 double const elastic_energy =
245 degradation * strain_energy_tensile + strain_energy_compressive;
246
247 KelvinVector sigma_real = degradation * sigma_tensile + sigma_compressive;
248 KelvinMatrix const D = degradation * C_tensile + C_compressive;
249
250 return std::make_tuple(sigma_real, sigma_tensile, D, strain_energy_tensile,
251 elastic_energy, C_tensile, C_compressive);
252}
253
254template <int DisplacementDim>
256 DisplacementDim> /* sigma_real */,
258 DisplacementDim> /* sigma_tensile */,
260 double /* strain_energy_tensile */, double /* elastic_energy */,
262 DisplacementDim> /* C_tensile */,
264 DisplacementDim> /* C_compressive
265 */
266 >
268 double const degradation,
269 double const bulk_modulus,
270 double const mu,
272{
273 static int const KelvinVectorSize =
275 using KelvinVector =
277 using KelvinMatrix =
280 // calculation of deviatoric parts
281 auto const& P_dev = Invariants::deviatoric_projection;
282 KelvinVector const epsd_curr = P_dev * eps;
283
284 // Hydrostatic part for the stress and the tangent.
285 double const eps_curr_trace = Invariants::trace(eps);
286
287 KelvinMatrix C_tensile = KelvinMatrix::Zero();
288 KelvinMatrix C_compressive = KelvinMatrix::Zero();
289
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 +
296 2 * mu * epsd_curr;
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;
300
301 KelvinMatrix const D = degradation * C_tensile + C_compressive;
302
303 return std::make_tuple(sigma_real, sigma_tensile, D, strain_energy_tensile,
304 elastic_energy, C_tensile, C_compressive);
305}
306
307template <int DisplacementDim>
309 DisplacementDim> /* sigma_real */,
311 DisplacementDim> /* sigma_tensile */,
313 double /* strain_energy_tensile */, double /* elastic_energy */,
315 DisplacementDim> /* C_tensile */,
317 DisplacementDim> /* C_compressive
318 */
319 >
321 double const degradation,
322 double const bulk_modulus,
323 double const mu,
325{
326 static int const KelvinVectorSize =
328 using KelvinVector =
330 using KelvinMatrix =
333 // calculation of deviatoric parts
334 auto const& P_dev = Invariants::deviatoric_projection;
335 KelvinVector const epsd_curr = P_dev * eps;
336
337 // Hydrostatic part for the stress and the tangent.
338 double const eps_curr_trace = Invariants::trace(eps);
339
340 KelvinMatrix C_tensile = KelvinMatrix::Zero();
341 KelvinMatrix C_compressive = KelvinMatrix::Zero();
342
343 KelvinVector const sigma_tensile =
344 bulk_modulus * eps_curr_trace * Invariants::identity2 +
345 2 * mu * epsd_curr;
346 // The maximum principal stress
347 auto sigma_tensor =
349 Eigen::SelfAdjointEigenSolver<decltype(sigma_tensor)>
350 selfAdjoint_eigen_solver(sigma_tensor);
351 double const max_principle_stress =
352 selfAdjoint_eigen_solver
353 .eigenvalues()[2]; // selfAdjoint_eigen_solver function gives an
354 // ascending sequence.
355 //
356 double const E0 = 9.0 * mu * bulk_modulus / (3.0 * bulk_modulus + mu);
357 auto hs = [&](double const v) { return heaviside(v); };
358 //
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);
365 //
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;
369
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);
373}
374
375} // namespace Phasefield
376} // namespace Solids
377} // 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