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