OGS
Lubby2.cpp
Go to the documentation of this file.
1
10#include "Lubby2.h"
11
12#include <Eigen/LU>
13
15
16namespace MPL = MaterialPropertyLib;
17
18namespace MaterialLib
19{
20namespace Solids
21{
22namespace Lubby2
23{
26template <int DisplacementDim>
27Eigen::Matrix<double, Lubby2<DisplacementDim>::JacobianResidualSize,
30{
31 Eigen::Matrix<double, Lubby2<DisplacementDim>::JacobianResidualSize,
33 dGdE =
34 Eigen::Matrix<double, Lubby2<DisplacementDim>::JacobianResidualSize,
36 dGdE.template topLeftCorner<Lubby2<DisplacementDim>::KelvinVectorSize,
38 .diagonal()
39 .setConstant(-2.);
40 return dGdE;
41}
42
43template <int DisplacementDim, typename LinearSolver>
45 double const GM0, double const KM0, LinearSolver const& linear_solver)
46{
47 // Calculate dGdE for time step
48 auto const dGdE = calculatedGdEBurgers<DisplacementDim>();
49
50 // Consistent tangent from local Newton iteration of material
51 // functionals.
52 // Only the upper left block is relevant for the global tangent.
53 static int const KelvinVectorSize =
55 using KelvinMatrix =
57
58 KelvinMatrix const dzdE =
59 linear_solver.solve(-dGdE)
60 .template topLeftCorner<KelvinVectorSize, KelvinVectorSize>();
61
63 auto const& P_sph = Invariants::spherical_projection;
64 auto const& P_dev = Invariants::deviatoric_projection;
65
66 KelvinMatrix C = GM0 * dzdE * P_dev + 3. * KM0 * P_sph;
67 return C;
68};
69
70template <int DisplacementDim>
71std::optional<std::tuple<typename Lubby2<DisplacementDim>::KelvinVector,
72 std::unique_ptr<typename MechanicsBase<
73 DisplacementDim>::MaterialStateVariables>,
76 MaterialPropertyLib::VariableArray const& variable_array_prev,
77 MaterialPropertyLib::VariableArray const& variable_array, double const t,
78 ParameterLib::SpatialPosition const& x, double const dt,
80 material_state_variables) const
81{
82 auto const& eps_m = std::get<MPL::SymmetricTensor<DisplacementDim>>(
83 variable_array.mechanical_strain);
84 auto const& eps_m_prev = std::get<MPL::SymmetricTensor<DisplacementDim>>(
85 variable_array_prev.mechanical_strain);
86 auto const& sigma_prev = std::get<MPL::SymmetricTensor<DisplacementDim>>(
87 variable_array_prev.stress);
88
90
91 assert(dynamic_cast<MaterialStateVariables const*>(
92 &material_state_variables) != nullptr);
94 static_cast<MaterialStateVariables const&>(material_state_variables));
96
97 auto local_lubby2_properties =
99
100 // calculation of deviatoric parts
101 auto const& P_dev = Invariants::deviatoric_projection;
102 KelvinVector const eps_m_d_i = P_dev * eps_m;
103 KelvinVector const eps_m_d_t = P_dev * eps_m_prev;
104
105 // initial guess as elastic predictor.
106 KelvinVector sigd_j = 2.0 * (eps_m_d_i - state.eps_M_t - state.eps_K_t);
107 // Note: sigd_t contains dimensionless stresses!
108 KelvinVector sigd_t = P_dev * sigma_prev / local_lubby2_properties.GM0;
109
110 // Calculate effective stress and update material properties
111 double sig_eff = Invariants::equivalentStress(sigd_j);
112 local_lubby2_properties.update(sig_eff);
113
114 using LocalJacobianMatrix =
115 Eigen::Matrix<double, KelvinVectorSize * 3, KelvinVectorSize * 3,
116 Eigen::RowMajor>;
117
118 // Linear solver for the newton loop is required after the loop with the
119 // same matrix. This saves one decomposition.
120 Eigen::FullPivLU<LocalJacobianMatrix> linear_solver;
121
122 // Different solvers are available for the solution of the local system.
123 // TODO Make the following choice of linear solvers available from the
124 // input file configuration:
125 // K_loc.partialPivLu().solve(-res_loc);
126 // K_loc.fullPivLu().solve(-res_loc);
127 // K_loc.householderQr().solve(-res_loc);
128 // K_loc.colPivHouseholderQr().solve(res_loc);
129 // K_loc.fullPivHouseholderQr().solve(-res_loc);
130 // K_loc.llt().solve(-res_loc);
131 // K_loc.ldlt().solve(-res_loc);
132
133 LocalJacobianMatrix K_loc;
134 { // Local Newton solver
135 using LocalResidualVector =
136 Eigen::Matrix<double, KelvinVectorSize * 3, 1>;
137
138 auto const update_residual = [&](LocalResidualVector& residual)
139 {
140 calculateResidualBurgers(dt, eps_m_d_i, eps_m_d_t, sigd_j, sigd_t,
141 state.eps_K_j, state.eps_K_t,
142 state.eps_M_j, state.eps_M_t, residual,
143 local_lubby2_properties);
144 };
145
146 auto const update_jacobian = [&](LocalJacobianMatrix& jacobian)
147 {
148 calculateJacobianBurgers(
149 t, x, dt, jacobian, sig_eff, sigd_j, state.eps_K_j,
150 local_lubby2_properties); // for solution dependent Jacobians
151 };
152
153 auto const update_solution = [&](LocalResidualVector const& increment)
154 {
155 // increment solution vectors
156 sigd_j.noalias() += increment.template segment<KelvinVectorSize>(
157 KelvinVectorSize * 0);
158 state.eps_K_j.noalias() +=
159 increment.template segment<KelvinVectorSize>(KelvinVectorSize *
160 1);
161 state.eps_M_j.noalias() +=
162 increment.template segment<KelvinVectorSize>(KelvinVectorSize *
163 2);
164
165 // Calculate effective stress and update material properties
167 KelvinVectorSize>::equivalentStress(sigd_j);
168 local_lubby2_properties.update(sig_eff);
169 };
170
171 auto newton_solver = NumLib::NewtonRaphson(
172 linear_solver, update_jacobian, update_residual, update_solution,
173 _nonlinear_solver_parameters);
174
175 auto const success_iterations = newton_solver.solve(K_loc);
176
177 if (!success_iterations)
178 {
179 return {};
180 }
181
182 // If the Newton loop didn't run, the linear solver will not be
183 // initialized.
184 // This happens usually for the first iteration of the first timestep.
185 if (*success_iterations == 0)
186 {
187 linear_solver.compute(K_loc);
188 }
189 }
190
191 KelvinMatrix C =
192 tangentStiffnessA<DisplacementDim>(local_lubby2_properties.GM0,
193 local_lubby2_properties.KM0,
194 linear_solver);
195
196 // Hydrostatic part for the stress and the tangent.
197 double const delta_eps_m_trace = Invariants::trace(eps_m - eps_m_prev);
198 double const sigma_trace_prev = Invariants::trace(sigma_prev);
199 KelvinVector const sigma =
200 local_lubby2_properties.GM0 * sigd_j +
201 (local_lubby2_properties.KM0 * delta_eps_m_trace +
202 sigma_trace_prev / 3.) *
203 Invariants::identity2;
204 return {std::make_tuple(
205 sigma,
206 std::unique_ptr<
208 new MaterialStateVariables{state}},
209 C)};
210}
211
212template <int DisplacementDim>
214 const double dt,
215 const KelvinVector& strain_curr,
216 const KelvinVector& strain_t,
217 const KelvinVector& stress_curr,
218 const KelvinVector& stress_t,
219 const KelvinVector& strain_Kel_curr,
220 const KelvinVector& strain_Kel_t,
221 const KelvinVector& strain_Max_curr,
222 const KelvinVector& strain_Max_t,
223 ResidualVector& res,
225{
226 // calculate stress residual
227 res.template segment<KelvinVectorSize>(0).noalias() =
228 (stress_curr - stress_t) -
229 2. * ((strain_curr - strain_t) - (strain_Kel_curr - strain_Kel_t) -
230 (strain_Max_curr - strain_Max_t));
231
232 // calculate Kelvin strain residual
233 res.template segment<KelvinVectorSize>(KelvinVectorSize).noalias() =
234 (strain_Kel_curr - strain_Kel_t) -
235 dt / (2. * properties.etaK) *
236 (properties.GM0 * stress_curr -
237 2. * properties.GK * strain_Kel_curr);
238
239 // calculate Maxwell strain residual
240 res.template segment<KelvinVectorSize>(2 * KelvinVectorSize).noalias() =
241 (strain_Max_curr - strain_Max_t) -
242 dt * 0.5 * properties.GM0 / properties.etaM * stress_curr;
243}
244
245template <int DisplacementDim>
247 double const t,
249 const double dt,
250 JacobianMatrix& Jac,
251 double s_eff,
252 const KelvinVector& sig_i,
253 const KelvinVector& eps_K_i,
255{
256 Jac.setZero();
257
258 // build G_11
259 Jac.template block<KelvinVectorSize, KelvinVectorSize>(0, 0)
260 .diagonal()
261 .setConstant(1.);
262
263 // build G_12
264 Jac.template block<KelvinVectorSize, KelvinVectorSize>(0, KelvinVectorSize)
265 .diagonal()
266 .setConstant(2.);
267
268 // build G_13
269 Jac.template block<KelvinVectorSize, KelvinVectorSize>(0,
270 2 * KelvinVectorSize)
271 .diagonal()
272 .setConstant(2.);
273
274 // build G_21
275 Jac.template block<KelvinVectorSize, KelvinVectorSize>(KelvinVectorSize, 0)
276 .noalias() =
277 -0.5 * dt * properties.GM0 / properties.etaK * KelvinMatrix::Identity();
278 if (s_eff > 0.)
279 {
280 KelvinVector const eps_K_aid =
281 1. / (properties.etaK * properties.etaK) *
282 (properties.GM0 * sig_i - 2. * properties.GK * eps_K_i);
283
284 KelvinVector const dG_K = 1.5 * _mp.mK(t, x)[0] * properties.GK *
285 properties.GM0 / s_eff * sig_i;
286 KelvinVector const dmu_vK = 1.5 * _mp.mvK(t, x)[0] * properties.GM0 *
287 properties.etaK / s_eff * sig_i;
288 Jac.template block<KelvinVectorSize, KelvinVectorSize>(KelvinVectorSize,
289 0)
290 .noalias() += 0.5 * dt * eps_K_aid * dmu_vK.transpose() +
291 dt / properties.etaK * eps_K_i * dG_K.transpose();
292 }
293
294 // build G_22
295 Jac.template block<KelvinVectorSize, KelvinVectorSize>(KelvinVectorSize,
296 KelvinVectorSize)
297 .diagonal()
298 .setConstant(1. + dt * properties.GK / properties.etaK);
299
300 // nothing to do for G_23
301
302 // build G_31
303 Jac.template block<KelvinVectorSize, KelvinVectorSize>(2 * KelvinVectorSize,
304 0)
305 .noalias() =
306 -0.5 * dt * properties.GM0 / properties.etaM * KelvinMatrix::Identity();
307 if (s_eff > 0.)
308 {
309 KelvinVector const dmu_vM = 1.5 * _mp.mvM(t, x)[0] * properties.GM0 *
310 properties.etaM / s_eff * sig_i;
311 Jac.template block<KelvinVectorSize, KelvinVectorSize>(
312 2 * KelvinVectorSize, 0)
313 .noalias() += 0.5 * dt * properties.GM0 /
314 (properties.etaM * properties.etaM) * sig_i *
315 dmu_vM.transpose();
316 }
317
318 // nothing to do for G_32
319
320 // build G_33
321 Jac.template block<KelvinVectorSize, KelvinVectorSize>(2 * KelvinVectorSize,
322 2 * KelvinVectorSize)
323 .diagonal()
324 .setConstant(1.);
325}
326
327template class Lubby2<2>;
328template class Lubby2<3>;
329
330} // namespace Lubby2
331} // namespace Solids
332} // namespace MaterialLib
std::optional< std::tuple< KelvinVector, std::unique_ptr< typename MechanicsBase< DisplacementDim >::MaterialStateVariables >, KelvinMatrix > > integrateStress(MaterialPropertyLib::VariableArray const &variable_array_prev, MaterialPropertyLib::VariableArray const &variable_array, double const t, ParameterLib::SpatialPosition const &x, double const dt, typename MechanicsBase< DisplacementDim >::MaterialStateVariables const &material_state_variables) const override
Definition Lubby2.cpp:75
MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > KelvinMatrix
Definition Lubby2.h:163
Eigen::Matrix< double, JacobianResidualSize, 1 > ResidualVector
Definition Lubby2.h:169
Eigen::Matrix< double, JacobianResidualSize, JacobianResidualSize, Eigen::RowMajor > JacobianMatrix
Definition Lubby2.h:170
static int const KelvinVectorSize
Definition Lubby2.h:159
MathLib::KelvinVector::KelvinVectorType< DisplacementDim > KelvinVector
Definition Lubby2.h:161
void calculateJacobianBurgers(double const t, ParameterLib::SpatialPosition const &x, double const dt, JacobianMatrix &Jac, double s_eff, const KelvinVector &sig_i, const KelvinVector &eps_K_i, detail::LocalLubby2Properties< DisplacementDim > const &properties) const
Calculates the 18x18 Jacobian.
Definition Lubby2.cpp:246
void calculateResidualBurgers(double const dt, const KelvinVector &strain_curr, const KelvinVector &strain_t, const KelvinVector &stress_curr, const KelvinVector &stress_t, const KelvinVector &strain_Kel_curr, const KelvinVector &strain_Kel_t, const KelvinVector &strain_Max_curr, const KelvinVector &strain_Max_t, ResidualVector &res, detail::LocalLubby2Properties< DisplacementDim > const &properties) const
Calculates the 18x1 residual vector.
Definition Lubby2.cpp:213
MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > tangentStiffnessA(double const GM0, double const KM0, LinearSolver const &linear_solver)
Definition Lubby2.cpp:44
Eigen::Matrix< double, Lubby2< DisplacementDim >::JacobianResidualSize, Lubby2< DisplacementDim >::KelvinVectorSize > calculatedGdEBurgers()
Definition Lubby2.cpp:29
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor > KelvinMatrixType