OGS
MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim > Class Template Referencefinal

Detailed Description

template<int DisplacementDim>
class MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >

Definition at line 98 of file Lubby2.h.

#include <Lubby2.h>

Inheritance diagram for MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >:
[legend]
Collaboration diagram for MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >:
[legend]

Classes

struct  MaterialStateVariables

Public Types

using KelvinVector
using KelvinMatrix
using ResidualVector = Eigen::Matrix<double, JacobianResidualSize, 1>
using JacobianMatrix
Public Types inherited from MaterialLib::Solids::MechanicsBase< DisplacementDim >
using KelvinVector
using KelvinMatrix

Public Member Functions

std::unique_ptr< typename MechanicsBase< DisplacementDim >::MaterialStateVariables > createMaterialStateVariables () const override
 Lubby2 (NumLib::NewtonRaphsonSolverParameters nonlinear_solver_parameters, Lubby2MaterialProperties &material_properties)
double computeFreeEnergyDensity (double const t, ParameterLib::SpatialPosition const &x, double const dt, KelvinVector const &eps, KelvinVector const &sigma, typename MechanicsBase< DisplacementDim >::MaterialStateVariables const &material_state_variables) const override
double getBulkModulus (double const t, ParameterLib::SpatialPosition const &x, KelvinMatrix const *const) const override
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
Public Member Functions inherited from MaterialLib::Solids::MechanicsBase< DisplacementDim >
virtual void initializeInternalStateVariables (double const, ParameterLib::SpatialPosition const &, typename MechanicsBase< DisplacementDim >::MaterialStateVariables &) const
virtual std::optional< std::tuple< KelvinVector, std::unique_ptr< MaterialStateVariables >, KelvinMatrix > > integrateStress (MaterialPropertyLib::VariableArray const &variable_array_prev, MaterialPropertyLib::VariableArray const &variable_array, double const t, ParameterLib::SpatialPosition const &x, double const dt, MaterialStateVariables const &material_state_variables) const =0
virtual std::vector< InternalVariablegetInternalVariables () const
virtual ConstitutiveModel getConstitutiveModel () const
 Gets the type of constitutive model.
virtual double getTemperatureRelatedCoefficient (double const, double const, ParameterLib::SpatialPosition const &, double const, double const) const
virtual double computeFreeEnergyDensity (double const t, ParameterLib::SpatialPosition const &x, double const dt, KelvinVector const &eps, KelvinVector const &sigma, MaterialStateVariables const &material_state_variables) const =0
virtual ~MechanicsBase ()=default

Static Public Attributes

static int const KelvinVectorSize
static int const JacobianResidualSize

Private Member Functions

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.
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.

Private Attributes

NumLib::NewtonRaphsonSolverParameters const _nonlinear_solver_parameters
Lubby2MaterialProperties _mp

Member Typedef Documentation

◆ JacobianMatrix

template<int DisplacementDim>
using MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >::JacobianMatrix
Initial value:
Eigen::Matrix<double,
Eigen::RowMajor>

Definition at line 164 of file Lubby2.h.

◆ KelvinMatrix

template<int DisplacementDim>
using MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >::KelvinMatrix
Initial value:
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor > KelvinMatrixType

Definition at line 157 of file Lubby2.h.

◆ KelvinVector

template<int DisplacementDim>
using MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >::KelvinVector
Initial value:
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType

Definition at line 155 of file Lubby2.h.

◆ ResidualVector

template<int DisplacementDim>
using MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >::ResidualVector = Eigen::Matrix<double, JacobianResidualSize, 1>

Definition at line 163 of file Lubby2.h.

Constructor & Destructor Documentation

◆ Lubby2()

template<int DisplacementDim>
MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >::Lubby2 ( NumLib::NewtonRaphsonSolverParameters nonlinear_solver_parameters,
Lubby2MaterialProperties & material_properties )
inlineexplicit

Definition at line 170 of file Lubby2.h.

References _mp, and _nonlinear_solver_parameters.

Member Function Documentation

◆ calculateJacobianBurgers()

template<int DisplacementDim>
void MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >::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
private

Calculates the 18x18 Jacobian.

Definition at line 240 of file Lubby2.cpp.

249{
250 Jac.setZero();
251
252 // build G_11
254 .diagonal()
255 .setConstant(1.);
256
257 // build G_12
259 .diagonal()
260 .setConstant(2.);
261
262 // build G_13
265 .diagonal()
266 .setConstant(2.);
267
268 // build G_21
270 .noalias() =
271 -0.5 * dt * properties.GM0 / properties.etaK * KelvinMatrix::Identity();
272 if (s_eff > 0.)
273 {
274 KelvinVector const eps_K_aid =
275 1. / (properties.etaK * properties.etaK) *
276 (properties.GM0 * sig_i - 2. * properties.GK * eps_K_i);
277
278 KelvinVector const dG_K = 1.5 * _mp.mK(t, x)[0] * properties.GK *
279 properties.GM0 / s_eff * sig_i;
280 KelvinVector const dmu_vK = 1.5 * _mp.mvK(t, x)[0] * properties.GM0 *
281 properties.etaK / s_eff * sig_i;
283 0)
284 .noalias() += 0.5 * dt * eps_K_aid * dmu_vK.transpose() +
285 dt / properties.etaK * eps_K_i * dG_K.transpose();
286 }
287
288 // build G_22
291 .diagonal()
292 .setConstant(1. + dt * properties.GK / properties.etaK);
293
294 // nothing to do for G_23
295
296 // build G_31
298 0)
299 .noalias() =
300 -0.5 * dt * properties.GM0 / properties.etaM * KelvinMatrix::Identity();
301 if (s_eff > 0.)
302 {
303 KelvinVector const dmu_vM = 1.5 * _mp.mvM(t, x)[0] * properties.GM0 *
304 properties.etaM / s_eff * sig_i;
306 2 * KelvinVectorSize, 0)
307 .noalias() += 0.5 * dt * properties.GM0 /
308 (properties.etaM * properties.etaM) * sig_i *
309 dmu_vM.transpose();
310 }
311
312 // nothing to do for G_32
313
314 // build G_33
317 .diagonal()
318 .setConstant(1.);
319}
static int const KelvinVectorSize
Definition Lubby2.h:153
MathLib::KelvinVector::KelvinVectorType< DisplacementDim > KelvinVector
Definition Lubby2.h:155

References _mp, MaterialLib::Solids::Lubby2::detail::LocalLubby2Properties< DisplacementDim >::etaK, MaterialLib::Solids::Lubby2::detail::LocalLubby2Properties< DisplacementDim >::etaM, MaterialLib::Solids::Lubby2::detail::LocalLubby2Properties< DisplacementDim >::GK, MaterialLib::Solids::Lubby2::detail::LocalLubby2Properties< DisplacementDim >::GM0, and KelvinVectorSize.

Referenced by integrateStress().

◆ calculateResidualBurgers()

template<int DisplacementDim>
void MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >::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
private

Calculates the 18x1 residual vector.

Definition at line 207 of file Lubby2.cpp.

219{
220 // calculate stress residual
221 res.template segment<KelvinVectorSize>(0).noalias() =
225
226 // calculate Kelvin strain residual
229 dt / (2. * properties.etaK) *
230 (properties.GM0 * stress_curr -
231 2. * properties.GK * strain_Kel_curr);
232
233 // calculate Maxwell strain residual
234 res.template segment<KelvinVectorSize>(2 * KelvinVectorSize).noalias() =
236 dt * 0.5 * properties.GM0 / properties.etaM * stress_curr;
237}

References MaterialLib::Solids::Lubby2::detail::LocalLubby2Properties< DisplacementDim >::etaK, MaterialLib::Solids::Lubby2::detail::LocalLubby2Properties< DisplacementDim >::etaM, MaterialLib::Solids::Lubby2::detail::LocalLubby2Properties< DisplacementDim >::GK, MaterialLib::Solids::Lubby2::detail::LocalLubby2Properties< DisplacementDim >::GM0, and KelvinVectorSize.

Referenced by integrateStress().

◆ computeFreeEnergyDensity()

template<int DisplacementDim>
double MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >::computeFreeEnergyDensity ( double const t,
ParameterLib::SpatialPosition const & x,
double const dt,
KelvinVector const & eps,
KelvinVector const & sigma,
typename MechanicsBase< DisplacementDim >::MaterialStateVariables const & material_state_variables ) const
inlineoverride

Definition at line 178 of file Lubby2.h.

186 {
187 assert(dynamic_cast<MaterialStateVariables const*>(
188 &material_state_variables) != nullptr);
191
192 auto const& eps_K = state.eps_K_j;
193 auto const& eps_K_prev = state.eps_K_t;
194
195 auto const& eps_M = state.eps_M_j;
196
199
200 // calculation of deviatoric parts
203 KelvinVector const epsd_i = P_dev * eps;
204
205 // initial guess as elastic predictor
206 KelvinVector sigd_j = 2.0 * (epsd_i - state.eps_M_t - state.eps_K_t);
207
208 // Calculate effective stress and update material properties
211
212 auto const& eta_K = local_lubby2_properties.etaK;
213
214 // This is specific to the backward Euler time scheme and needs to be
215 // updated if other time schemes are used.
216 return (eps - eps_K - eps_M).dot(sigma) / 2 +
217 eps_K.dot(sigma - eta_K * (eps_K - eps_K_prev) / dt) / 2;
218 }

References _mp, MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >::MaterialStateVariables::eps_K_j, MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >::MaterialStateVariables::eps_K_t, MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >::MaterialStateVariables::eps_M_j, and MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >::MaterialStateVariables::eps_M_t.

◆ createMaterialStateVariables()

template<int DisplacementDim>
std::unique_ptr< typename MechanicsBase< DisplacementDim >::MaterialStateVariables > MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >::createMaterialStateVariables ( ) const
inlineoverridevirtual

Polymorphic creator for MaterialStateVariables objects specific for a material model.

Reimplemented from MaterialLib::Solids::MechanicsBase< DisplacementDim >.

Definition at line 145 of file Lubby2.h.

◆ getBulkModulus()

template<int DisplacementDim>
double MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >::getBulkModulus ( double const t,
ParameterLib::SpatialPosition const & x,
KelvinMatrix const * const  ) const
inlineoverridevirtual

Reimplemented from MaterialLib::Solids::MechanicsBase< DisplacementDim >.

Definition at line 220 of file Lubby2.h.

223 {
224 return _mp.KM0(t, x)[0];
225 }

References _mp.

◆ integrateStress()

template<int DisplacementDim>
std::optional< std::tuple< typename Lubby2< DisplacementDim >::KelvinVector, std::unique_ptr< typename MechanicsBase< DisplacementDim >::MaterialStateVariables >, typename Lubby2< DisplacementDim >::KelvinMatrix > > MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >::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 at line 69 of file Lubby2.cpp.

75{
77 variable_array.mechanical_strain);
79 variable_array_prev.mechanical_strain);
81 variable_array_prev.stress);
82
84
85 assert(dynamic_cast<MaterialStateVariables const*>(
86 &material_state_variables) != nullptr);
89 state.setInitialConditions();
90
93
94 // calculation of deviatoric parts
98
99 // initial guess as elastic predictor.
100 KelvinVector sigd_j = 2.0 * (eps_m_d_i - state.eps_M_t - state.eps_K_t);
101 // Note: sigd_t contains dimensionless stresses!
103
104 // Calculate effective stress and update material properties
107
108 using LocalJacobianMatrix =
111
112 // Linear solver for the newton loop is required after the loop with the
113 // same matrix. This saves one decomposition.
115
116 // Different solvers are available for the solution of the local system.
117 // TODO Make the following choice of linear solvers available from the
118 // input file configuration:
119 // K_loc.partialPivLu().solve(-res_loc);
120 // K_loc.fullPivLu().solve(-res_loc);
121 // K_loc.householderQr().solve(-res_loc);
122 // K_loc.colPivHouseholderQr().solve(res_loc);
123 // K_loc.fullPivHouseholderQr().solve(-res_loc);
124 // K_loc.llt().solve(-res_loc);
125 // K_loc.ldlt().solve(-res_loc);
126
128 { // Local Newton solver
129 using LocalResidualVector =
131
133 {
135 state.eps_K_j, state.eps_K_t,
136 state.eps_M_j, state.eps_M_t, residual,
138 };
139
141 {
143 t, x, dt, jacobian, sig_eff, sigd_j, state.eps_K_j,
144 local_lubby2_properties); // for solution dependent Jacobians
145 };
146
147 auto const update_solution = [&](LocalResidualVector const& increment)
148 {
149 // increment solution vectors
150 sigd_j.noalias() += increment.template segment<KelvinVectorSize>(
151 KelvinVectorSize * 0);
152 state.eps_K_j.noalias() +=
154 1);
155 state.eps_M_j.noalias() +=
157 2);
158
159 // Calculate effective stress and update material properties
163 };
164
168
169 auto const success_iterations = newton_solver.solve(K_loc);
170
172 {
173 return {};
174 }
175
176 // If the Newton loop didn't run, the linear solver will not be
177 // initialized.
178 // This happens usually for the first iteration of the first timestep.
179 if (*success_iterations == 0)
180 {
181 linear_solver.compute(K_loc);
182 }
183 }
184
189
190 // Hydrostatic part for the stress and the tangent.
193 KelvinVector const sigma =
196 sigma_trace_prev / 3.) *
198 return {std::make_tuple(
199 sigma,
203 C)};
204}
MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > KelvinMatrix
Definition Lubby2.h:157
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:240
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:207
MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > tangentStiffnessA(double const GM0, double const KM0, LinearSolver const &linear_solver)
Definition Lubby2.cpp:38

References _mp, _nonlinear_solver_parameters, calculateJacobianBurgers(), calculateResidualBurgers(), MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >::MaterialStateVariables::eps_K_j, MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >::MaterialStateVariables::eps_K_t, MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >::MaterialStateVariables::eps_M_j, MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >::MaterialStateVariables::eps_M_t, KelvinVectorSize, MaterialPropertyLib::VariableArray::mechanical_strain, MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >::MaterialStateVariables::setInitialConditions(), MaterialPropertyLib::VariableArray::stress, and MaterialLib::Solids::Lubby2::tangentStiffnessA().

Member Data Documentation

◆ _mp

template<int DisplacementDim>
Lubby2MaterialProperties MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >::_mp
private

◆ _nonlinear_solver_parameters

template<int DisplacementDim>
NumLib::NewtonRaphsonSolverParameters const MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >::_nonlinear_solver_parameters
private

Definition at line 265 of file Lubby2.h.

Referenced by Lubby2(), and integrateStress().

◆ JacobianResidualSize

template<int DisplacementDim>
int const MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >::JacobianResidualSize
static
Initial value:

Definition at line 160 of file Lubby2.h.

◆ KelvinVectorSize

template<int DisplacementDim>
int const MaterialLib::Solids::Lubby2::Lubby2< DisplacementDim >::KelvinVectorSize
static

The documentation for this class was generated from the following files: