OGS
MaterialLib::Fracture::Coulomb::Coulomb< DisplacementDim > Class Template Referencefinal

Detailed Description

template<int DisplacementDim>
class MaterialLib::Fracture::Coulomb::Coulomb< DisplacementDim >

Definition at line 39 of file Coulomb.h.

#include <Coulomb.h>

Inheritance diagram for MaterialLib::Fracture::Coulomb::Coulomb< DisplacementDim >:
[legend]
Collaboration diagram for MaterialLib::Fracture::Coulomb::Coulomb< DisplacementDim >:
[legend]

Classes

struct  MaterialProperties
 Variables specific to the material model. More...

Public Member Functions

std::unique_ptr< typename FractureModelBase< DisplacementDim >::MaterialStateVariables > createMaterialStateVariables () override
 Coulomb (NumLib::NewtonRaphsonSolverParameters nonlinear_solver_parameters, double const penalty_aperture_cutoff, bool const tension_cutoff, MaterialProperties material_properties)
void computeConstitutiveRelation (double const t, ParameterLib::SpatialPosition const &x, double const aperture0, Eigen::Ref< Eigen::VectorXd const > sigma0, Eigen::Ref< Eigen::VectorXd const > w_prev, Eigen::Ref< Eigen::VectorXd const > w, Eigen::Ref< Eigen::VectorXd const > sigma_prev, Eigen::Ref< Eigen::VectorXd > sigma, Eigen::Ref< Eigen::MatrixXd > Kep, typename FractureModelBase< DisplacementDim >::MaterialStateVariables &material_state_variables) override
Public Member Functions inherited from MaterialLib::Fracture::FractureModelBase< DisplacementDim >
virtual ~FractureModelBase ()=default
virtual void computeConstitutiveRelation (double const t, ParameterLib::SpatialPosition const &x, double const aperture0, Eigen::Ref< Eigen::VectorXd const > sigma0, Eigen::Ref< Eigen::VectorXd const > w_prev, Eigen::Ref< Eigen::VectorXd const > w, Eigen::Ref< Eigen::VectorXd const > sigma_prev, Eigen::Ref< Eigen::VectorXd > sigma, Eigen::Ref< Eigen::MatrixXd > C, MaterialStateVariables &material_state_variables)=0

Private Attributes

NumLib::NewtonRaphsonSolverParameters const _nonlinear_solver_parameters
double const _penalty_aperture_cutoff
bool const _tension_cutoff
MaterialProperties _mp

Constructor & Destructor Documentation

◆ Coulomb()

template<int DisplacementDim>
MaterialLib::Fracture::Coulomb::Coulomb< DisplacementDim >::Coulomb ( NumLib::NewtonRaphsonSolverParameters nonlinear_solver_parameters,
double const penalty_aperture_cutoff,
bool const tension_cutoff,
MaterialProperties material_properties )
inlineexplicit

Member Function Documentation

◆ computeConstitutiveRelation()

template<int DisplacementDim>
void MaterialLib::Fracture::Coulomb::Coulomb< DisplacementDim >::computeConstitutiveRelation ( double const t,
ParameterLib::SpatialPosition const & x,
double const aperture0,
Eigen::Ref< Eigen::VectorXd const > sigma0,
Eigen::Ref< Eigen::VectorXd const > w_prev,
Eigen::Ref< Eigen::VectorXd const > w,
Eigen::Ref< Eigen::VectorXd const > sigma_prev,
Eigen::Ref< Eigen::VectorXd > sigma,
Eigen::Ref< Eigen::MatrixXd > Kep,
typename FractureModelBase< DisplacementDim >::MaterialStateVariables & material_state_variables )
override

Computation of the constitutive relation for the Mohr-Coulomb model.

Parameters
tcurrent time
xcurrent position in space
aperture0initial fracture's aperture
sigma0initial stress
w_prevfracture displacement at previous time step
wfracture displacement at current time step
sigma_prevstress at previous time step
sigmastress at current time step
Keptangent matrix for stress and fracture displacements
material_state_variablesmaterial state variables

Definition at line 44 of file Coulomb.cpp.

62{
64 &material_state_variables) != nullptr);
65
68
70
71 const int index_ns = DisplacementDim - 1;
72 double const aperture = w[index_ns] + aperture0;
73
75 { // Elastic tangent stiffness
77 for (int i = 0; i < index_ns; i++)
78 {
79 Ke(i, i) = mat.Ks;
80 }
81
82 Ke(index_ns, index_ns) = mat.Kn;
83 }
84
85 // Total plastic aperture compression
86 // NOTE: Initial condition sigma0 seems to be associated with an initial
87 // condition of the w0 = 0. Therefore the initial state is not associated
88 // with a plastic aperture change.
89 { // Exact elastic predictor
90 sigma.noalias() = Ke * (w - w_prev);
91
92 sigma.coeffRef(index_ns) *=
94 sigma.noalias() += sigma_prev;
95 }
96
97 // correction for an opening fracture
98 if (_tension_cutoff && sigma[DisplacementDim - 1] >= 0)
99 {
100 Kep.setZero();
101 sigma.setZero();
102 state.w_p = w;
103 material_state_variables.setTensileStress(true);
104 return;
105 }
106
107 auto yield_function = [&mat](Eigen::VectorXd const& s)
108 {
109 double const sigma_n = s[DisplacementDim - 1];
110 Eigen::VectorXd const sigma_s = s.head(DisplacementDim - 1);
111 double const mag_tau = sigma_s.norm(); // magnitude
112 return mag_tau + sigma_n * std::tan(mat.phi) - mat.c;
113 };
114
115 { // Exit if still in elastic range by checking the shear yield function.
116 double const Fs = yield_function(sigma);
117 material_state_variables.setShearYieldFunctionValue(Fs);
118 if (Fs < .0)
119 {
120 Kep = Ke;
123 return;
124 }
125 }
126
128 {
130 dFs_dS.template head<DisplacementDim - 1>().noalias() =
131 s.template head<DisplacementDim - 1>().normalized();
132 dFs_dS.coeffRef(DisplacementDim - 1) = std::tan(mat.phi);
133 return dFs_dS;
134 };
135
136 // plastic potential function: Qs = |tau| + Sn * tan da
138 {
140 dQs_dS.template head<DisplacementDim - 1>().noalias() =
141 s.template head<DisplacementDim - 1>().normalized();
142 dQs_dS.coeffRef(DisplacementDim - 1) = std::tan(mat.psi);
143 return dQs_dS;
144 };
145
146 { // Newton
147
152
155 solution << 0;
156
158 { residual[0] = yield_function(sigma); };
159
160 auto const update_jacobian = [&](JacobianMatrix& jacobian)
161 {
162 jacobian(0, 0) = -yield_function_derivative(sigma).transpose() *
164 };
165
166 auto const update_solution = [&](ResidualVectorType const& increment)
167 {
169 /*DBUG("analytical = {:g}",
170 Fs / (mat.Ks + mat.Kn * std::tan(mat.psi) * std::tan(mat.phi)))
171 */
172 state.w_p = state.w_p_prev +
174
175 sigma.noalias() = Ke * (w - w_prev - state.w_p + state.w_p_prev);
176
179 sigma.noalias() += sigma_prev;
180 };
181
185
186 auto const success_iterations = newton_solver.solve(jacobian);
187
189 {
191 "FractureModel/Coulomb local nonlinear solver didn't "
192 "converge.");
193 }
194
195 // Solution containing lambda is not needed; w_p and sigma already
196 // up to date.
197 }
198
199 { // Update material state shear yield function value.
200 double const Fs = yield_function(sigma);
201 material_state_variables.setShearYieldFunctionValue(Fs);
202 }
203
204 Ke(index_ns, index_ns) *=
207 Ke /
208 (yield_function_derivative(sigma).transpose() *
211}
double logPenaltyDerivative(double const aperture0, double const aperture, double const aperture_cutoff)
Definition LogPenalty.h:17

References _mp, _nonlinear_solver_parameters, _penalty_aperture_cutoff, _tension_cutoff, MaterialLib::Fracture::logPenaltyDerivative(), MaterialLib::Fracture::FractureModelBase< DisplacementDim >::MaterialStateVariables::setShearYieldFunctionValue(), MaterialLib::Fracture::FractureModelBase< DisplacementDim >::MaterialStateVariables::setTensileStress(), MaterialLib::Fracture::Coulomb::StateVariables< DisplacementDim >::w_p, and MaterialLib::Fracture::Coulomb::StateVariables< DisplacementDim >::w_p_prev.

◆ createMaterialStateVariables()

template<int DisplacementDim>
std::unique_ptr< typename FractureModelBase< DisplacementDim >::MaterialStateVariables > MaterialLib::Fracture::Coulomb::Coulomb< DisplacementDim >::createMaterialStateVariables ( )
inlineoverridevirtual

Polymorphic creator for MaterialStateVariables objects specific for a material model.

Implements MaterialLib::Fracture::FractureModelBase< DisplacementDim >.

Definition at line 44 of file Coulomb.h.

Member Data Documentation

◆ _mp

template<int DisplacementDim>
MaterialProperties MaterialLib::Fracture::Coulomb::Coulomb< DisplacementDim >::_mp
private

Definition at line 141 of file Coulomb.h.

Referenced by Coulomb(), and computeConstitutiveRelation().

◆ _nonlinear_solver_parameters

template<int DisplacementDim>
NumLib::NewtonRaphsonSolverParameters const MaterialLib::Fracture::Coulomb::Coulomb< DisplacementDim >::_nonlinear_solver_parameters
private

Definition at line 129 of file Coulomb.h.

Referenced by Coulomb(), and computeConstitutiveRelation().

◆ _penalty_aperture_cutoff

template<int DisplacementDim>
double const MaterialLib::Fracture::Coulomb::Coulomb< DisplacementDim >::_penalty_aperture_cutoff
private

Compressive normal displacements above this value will not enter the computation of the normal stiffness modulus of the fracture.

Note
Setting this to the initial aperture value allows negative apertures.

Definition at line 135 of file Coulomb.h.

Referenced by Coulomb(), and computeConstitutiveRelation().

◆ _tension_cutoff

template<int DisplacementDim>
bool const MaterialLib::Fracture::Coulomb::Coulomb< DisplacementDim >::_tension_cutoff
private

If set no resistance to open the fracture over the initial aperture is opposed.

Definition at line 139 of file Coulomb.h.

Referenced by Coulomb(), and computeConstitutiveRelation().


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