OGS
PhaseFieldBase.h
Go to the documentation of this file.
1
9
10#pragma once
11
12#include <Eigen/Core>
13#include <numbers>
14#include <optional>
15
16#include "BaseLib/Logging.h"
21
22namespace MaterialLib
23{
24namespace Solids
25{
26namespace Phasefield
27{
34
35template <int DisplacementDim>
37 std::string const& phasefield_model)
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};
56
62
63template <int DisplacementDim>
65 std::optional<std::string> const& softening_curve)
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};
84
94
95template <int DisplacementDim>
97 std::string const& energy_split_model)
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};
128
130{
131public:
132 virtual ~DegradationDerivative() = default;
133 virtual double degradation(double const d_ip,
134 double const k,
135 double const ls) = 0; /* degradation_df0 */
136 virtual double degradationDf1(double const d_ip, double const k,
137 double const ls) = 0;
138 virtual double degradationDf2(double const d_ip, double const k,
139 double const ls) = 0;
140};
141
143{
144public:
145 double degradation(double const d_ip,
146 double const k,
147 double const /* ls */) override
148 {
149 return d_ip * d_ip * (1. - k) + k;
150 };
151 double degradationDf1(double const d_ip, double const k,
152 double const /* ls */) override
153 {
154 return 2. * (1. - k) * d_ip;
155 };
156 double degradationDf2(double const /* d_ip */, double const k,
157 double const /* ls */) override
158 {
159 return 2. * (1. - k);
160 };
161};
162
164{
165private:
166 double const lch;
168
169public:
170 COHESIVE_DegradationDerivative(double const Parameter_lch,
171 SofteningCurve Parameter_softening_curve)
172 : lch(Parameter_lch), softening_curve(Parameter_softening_curve){};
173 double degradation(double const d_ip,
174 double const k,
175 double const ls) override
176 {
177 double const m1 = 4.0 * lch / acos(-1.) / ls;
178 switch (softening_curve)
179 {
181 {
182 double const m2 = std::pow(2., 5. / 3.) - 3.;
183 double const n = 2.5;
184 return std::pow(d_ip, n) /
185 (std::pow(d_ip, n) +
186 m1 * (1. - d_ip) * (1. + m2 * (1. - d_ip))) *
187 (1.0 - k) +
188 k;
189 }
190 default:
191 {
192 double const m2 = -0.5;
193 double const n = 2.;
194 return std::pow(d_ip, n) /
195 (std::pow(d_ip, n) +
196 m1 * (1. - d_ip) * (1. + m2 * (1. - d_ip))) *
197 (1. - k) +
198 k;
199 }
200 }
201 };
202 double degradationDf1(double const d_ip, double const k,
203 double const ls) override
204 {
205 double const m1 = 4.0 * lch / acos(-1.) / ls;
206 switch (softening_curve)
207 {
209 {
210 double const m2 = std::pow(2., 5. / 3.) - 3.;
211 double const n = 2.5;
212 double const a1 = std::pow(d_ip, n) +
213 m1 * (1. - d_ip) * (1. + m2 * (1. - d_ip));
214 double const a2 = n * std::pow(d_ip, n - 1.) -
215 2. * m1 * m2 * (1. - d_ip) - m1;
216 return (1. - k) *
217 (n * std::pow(d_ip, n - 1.) * a1 -
218 std::pow(d_ip, n) * a2) /
219 (a1 * a1);
220 }
221 default:
222 {
223 double const m2 = -0.5;
224 double const n = 2.;
225 double const a1 = std::pow(d_ip, n) +
226 m1 * (1. - d_ip) * (1. + m2 * (1. - d_ip));
227 double const a2 = n * std::pow(d_ip, n - 1.) -
228 2. * m1 * m2 * (1. - d_ip) - m1;
229 return (1. - k) *
230 (n * std::pow(d_ip, n - 1.) * a1 -
231 std::pow(d_ip, n) * a2) /
232 (a1 * a1);
233 }
234 }
235 };
236 double degradationDf2(double const d_ip, double const k,
237 double const ls) override
238 {
239 double const m1 = 4.0 * lch / acos(-1.) / ls;
240 switch (softening_curve)
241 {
243 {
244 double const m2 = std::pow(2., 5. / 3.) - 3.;
245 double const n = 2.5;
246 double a1 = std::pow(d_ip, n) +
247 m1 * (1. - d_ip) * (1. + m2 * (1. - d_ip));
248 double a2 = n * std::pow(d_ip, n - 1.) -
249 2. * m1 * m2 * (1. - d_ip) - m1;
250 return (1. - k) *
251 (2. * a2 * a2 * std::pow(d_ip, n) -
252 a1 * std::pow(d_ip, n) *
253 (2. * m1 * m2 +
254 n * std::pow(d_ip, n - 2.) * (n - 1.)) -
255 2. * a1 * a2 * n * std::pow(d_ip, n - 1.) +
256 a1 * a1 * n * std::pow(d_ip, n - 2.) * (n - 1.)) /
257 (std::pow(a1, 3));
258 }
259 default:
260 {
261 double const m2 = -0.5;
262 double const n = 2.;
263 double const a1 = std::pow(d_ip, n) +
264 m1 * (1. - d_ip) * (1. + m2 * (1. - d_ip));
265 double const a2 = n * std::pow(d_ip, n - 1.) -
266 2. * m1 * m2 * (1. - d_ip) - m1;
267 return (1. - k) *
268 (2. * a2 * a2 * std::pow(d_ip, n) -
269 a1 * std::pow(d_ip, n) *
270 (2. * m1 * m2 +
271 n * std::pow(d_ip, n - 2.) * (n - 1.)) -
272 2. * a1 * a2 * n * std::pow(d_ip, n - 1.) +
273 a1 * a1 * n * std::pow(d_ip, n - 2.) * (n - 1.)) /
274 (std::pow(a1, 3));
275 }
276 }
277 };
278};
279
280template <int DisplacementDim>
281std::unique_ptr<DegradationDerivative> creatDegradationDerivative(
282 PhaseFieldModel const& phasefield_model, double const lch,
283 SofteningCurve const& softening_curve)
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};
294
295template <typename T_DNDX, typename T_N, typename T_W, typename T_D,
296 typename T_LOCAL_JAC, typename T_LOCAL_RHS>
297void calculateCrackLocalJacobianAndResidual(T_DNDX& dNdx, T_N& N, T_W& w,
298 T_D& d, T_LOCAL_JAC& local_Jac,
299 T_LOCAL_RHS local_rhs,
300 double const gc, double const ls,
301 PhaseFieldModel& phasefield_model)
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};
345
346template <typename T_VECTOR, typename T_MATRIX, int DisplacementDim>
348 T_VECTOR& sigma, T_VECTOR& sigma_tensile, T_VECTOR& sigma_compressive,
349 T_VECTOR& eps_tensile, T_MATRIX& D, T_MATRIX& C_tensile,
350 T_MATRIX& C_compressive, double& strain_energy_tensile,
351 double& elastic_energy, double const degradation, T_VECTOR const& eps,
352 EnergySplitModel const& energy_split_model, double const t,
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.;
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.;
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};
435} // namespace Phasefield
436} // namespace Solids
437} // namespace MaterialLib
#define OGS_FATAL(...)
Definition Error.h:26
double degradation(double const d_ip, double const k, double const) override
double degradationDf1(double const d_ip, double const k, double const) override
double degradationDf2(double const, double const k, double const) override
double degradationDf1(double const d_ip, double const k, double const ls) override
double degradationDf2(double const d_ip, double const k, double const ls) override
double degradation(double const d_ip, double const k, double const ls) override
COHESIVE_DegradationDerivative(double const Parameter_lch, SofteningCurve Parameter_softening_curve)
virtual double degradationDf2(double const d_ip, double const k, double const ls)=0
virtual double degradation(double const d_ip, double const k, double const ls)=0
virtual double degradationDf1(double const d_ip, double const k, double const ls)=0
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)
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)
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)
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)
std::unique_ptr< DegradationDerivative > creatDegradationDerivative(PhaseFieldModel const &phasefield_model, double const lch, SofteningCurve const &softening_curve)
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)
PhaseFieldModel convertStringToPhaseFieldModel(std::string const &phasefield_model)
SofteningCurve convertStringToSofteningCurve(std::optional< std::string > const &softening_curve)
EnergySplitModel convertStringToEnergySplitModel(std::string const &energy_split_model)
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)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor > KelvinMatrixType