OGS
LinearElasticOrthotropicPhaseField.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
7
8namespace MaterialLib
9{
10namespace Solids
11{
12namespace Phasefield
13{
14
15template <int DisplacementDim>
17 DisplacementDim> /* eps_tensile */,
19 DisplacementDim> /* sigma_real */,
21 DisplacementDim> /* sigma_tensile */,
23 DisplacementDim> /* sigma_compressive */,
25 double /* strain_energy_tensile */, double /* elastic_energy */,
27 DisplacementDim> /* C_tensile */,
29 DisplacementDim> /* C_compressive
30 */
31 >
33 double const degradation,
36{
37 static constexpr int KelvinVectorSize =
39 using KelvinVector =
41 using KelvinMatrix =
44 // calculation of deviatoric parts
45 auto const& P_dev = Invariants::deviatoric_projection;
46
47 KelvinMatrix const C = C_ortho;
48 // calculation of square root of elasticity tensor
49 Eigen::SelfAdjointEigenSolver<KelvinMatrix> es(C);
50 KelvinMatrix const sqrt_C = es.operatorSqrt();
51 Eigen::SelfAdjointEigenSolver<KelvinMatrix> es_inverse(C.inverse());
52 KelvinMatrix const sqrt_C_inverse = es_inverse.operatorSqrt();
53
54 // strain in a transformed space
55 KelvinVector const epst = sqrt_C * eps;
56 double const epst_curr_trace = Invariants::trace(epst);
57
58 // projection tensors in transformed space
59 KelvinMatrix teps_p = KelvinMatrix::Zero();
60 KelvinMatrix teps_n = KelvinMatrix::Zero();
61 if (epst_curr_trace >= 0) /* QQQ */
62 {
63 teps_p.template topLeftCorner<3, 3>().setConstant(1. / 3);
64 }
65 else
66 {
67 teps_n.template topLeftCorner<3, 3>().setConstant(1. / 3);
68 }
69
70 teps_p.noalias() += P_dev * KelvinMatrix::Identity();
71
72 // strain tensile / compressive
73 KelvinVector const eps_tensile = sqrt_C_inverse * (teps_p * sqrt_C) * eps;
74 KelvinVector const eps_compressive =
75 sqrt_C_inverse * (teps_n * sqrt_C) * eps;
76
77 // projection tensors in original space
78 KelvinMatrix const der_eps_p = sqrt_C_inverse * (teps_p * sqrt_C);
79 KelvinMatrix const der_eps_n = sqrt_C_inverse * (teps_n * sqrt_C);
80 // C tensile / compressive
81 KelvinMatrix const C_tensile = der_eps_p.transpose() * C * der_eps_p;
82 KelvinMatrix const C_compressive = der_eps_n.transpose() * C * der_eps_n;
83
84 // stress tensile / compressive
85 KelvinVector const sigma_tensile = C_tensile * eps;
86 KelvinVector const sigma_compressive = C_compressive * eps;
87
88 // decomposition of strain energy density
89 double const strain_energy_tensile =
90 0.5 * sigma_tensile.adjoint() * eps_tensile;
91 double const strain_energy_compressive =
92 0.5 * sigma_compressive.adjoint() * eps_compressive;
93
94 double const elastic_energy =
95 degradation * strain_energy_tensile + strain_energy_compressive;
96
97 KelvinVector const sigma_real =
98 degradation * sigma_tensile + sigma_compressive;
99 KelvinMatrix const D = degradation * C_tensile + C_compressive;
100
101 return std::make_tuple(eps_tensile, sigma_real, sigma_tensile,
102 sigma_compressive, D, strain_energy_tensile,
103 elastic_energy, C_tensile, C_compressive);
104}
105
106template <int DisplacementDim>
108 DisplacementDim> /* eps_tensile */,
110 DisplacementDim> /* sigma_real */,
112 DisplacementDim> /* sigma_tensile */,
114 double /* strain_energy_tensile */, double /* elastic_energy */,
116 DisplacementDim> /* C_tensile */,
118 DisplacementDim> /* C_compressive
119 */
120 >
122 double const degradation,
125{
126 using KelvinVector =
128 using KelvinMatrix =
130
131 KelvinMatrix const C = C_ortho;
132
133 // calculation of square root of elasticity tensor
134 Eigen::SelfAdjointEigenSolver<KelvinMatrix> es(C);
135 KelvinMatrix const sqrt_C = es.operatorSqrt();
136 Eigen::SelfAdjointEigenSolver<KelvinMatrix> es_inverse(C.inverse());
137 KelvinMatrix const sqrt_C_inverse = es_inverse.operatorSqrt();
138
139 // strain in a transformed space
140 KelvinVector const epst = sqrt_C * eps;
141
142 // strain tensor in 3D
143 Eigen::Matrix3d const eps_3D =
145
146 Eigen::SelfAdjointEigenSolver<Eigen::Matrix3d> const es_eps_3D(eps_3D);
147 Eigen::Vector3d const eigen_values_eps_3D = es_eps_3D.eigenvalues().real();
148 Eigen::Matrix3d const eigen_vectors_eps_3D = es_eps_3D.eigenvectors();
149 Eigen::Matrix3d const E1_eigenp =
150 eigen_vectors_eps_3D.col(0) * eigen_vectors_eps_3D.col(0).transpose();
151
152 Eigen::Matrix3d const E3_eigenp =
153 eigen_vectors_eps_3D.col(2) * eigen_vectors_eps_3D.col(2).transpose();
154
155 KelvinVector const E1_vec =
157
158 KelvinVector const E3_vec =
160
161 KelvinMatrix const E1oE1 = E1_vec * E1_vec.transpose();
162 KelvinMatrix const E3oE3 = E3_vec * E3_vec.transpose();
163
164 KelvinMatrix I_p = KelvinMatrix::Zero();
165 KelvinMatrix I_n = KelvinMatrix::Zero();
166
167 KelvinMatrix const I_S = KelvinMatrix::Identity();
168 if (DisplacementDim == 2)
169 {
170 if (std::abs(eigen_values_eps_3D(2) - eigen_values_eps_3D(0)) >
171 std::numeric_limits<double>::epsilon())
172 {
173 I_p = (macaulayTensile(eigen_values_eps_3D(0)) -
174 macaulayTensile(eigen_values_eps_3D(2))) /
175 (eigen_values_eps_3D(0) - eigen_values_eps_3D(2)) *
176 (I_S - (E1oE1 + E3oE3)) +
177 (heaviside(eigen_values_eps_3D(0)) * E1oE1 +
178 heaviside(eigen_values_eps_3D(2)) * E3oE3);
179 I_n = (macaulayCompressive(eigen_values_eps_3D(0)) -
180 macaulayCompressive(eigen_values_eps_3D(2))) /
181 (eigen_values_eps_3D(0) - eigen_values_eps_3D(2)) *
182 (I_S - (E1oE1 + E3oE3)) +
183 (heaviside(-eigen_values_eps_3D(0)) * E1oE1 +
184 heaviside(-eigen_values_eps_3D(2)) * E3oE3);
185 }
186 else
187 {
188 I_p = heaviside(eigen_values_eps_3D(0)) * I_S;
189 I_n = heaviside(-eigen_values_eps_3D(0)) * I_S;
190 }
191 }
192
193 // strain tensile / compressive
194 KelvinVector const eps_tensile = sqrt_C_inverse * (I_p * sqrt_C) * eps;
195 KelvinVector const eps_compressive = sqrt_C_inverse * (I_n * sqrt_C) * eps;
196
197 // C tensile / compressive
198 KelvinMatrix const C_tensile = sqrt_C * I_p * sqrt_C;
199 KelvinMatrix const C_compressive = sqrt_C * I_n * sqrt_C;
200
201 // stress tensile / compressive
202 KelvinVector const sigma_tensile = C * eps_tensile;
203 KelvinVector const sigma_compressive = C * eps_compressive;
204
205 // decomposition of strain energy density
206 double const strain_energy_tensile =
207 0.5 * sigma_tensile.adjoint() * eps_tensile;
208 double const strain_energy_compressive =
209 0.5 * sigma_compressive.adjoint() * eps_compressive;
210
211 double const elastic_energy =
212 degradation * strain_energy_tensile + strain_energy_compressive;
213
214 KelvinVector const sigma_real =
215 degradation * sigma_tensile + sigma_compressive;
216
217 KelvinMatrix const D = degradation * C_tensile + C_compressive;
218
219 return std::make_tuple(eps_tensile, sigma_real, sigma_tensile, D,
220 strain_energy_tensile, elastic_energy, C_tensile,
221 C_compressive);
222}
223
224} // namespace Phasefield
225} // namespace Solids
226} // namespace MaterialLib
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)
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)
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