OGS
Coulomb.cpp
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#include "Coulomb.h"
5
6#include <Eigen/LU>
7#include <boost/math/constants/constants.hpp>
8
9#include "BaseLib/Error.h"
10#include "LogPenalty.h"
11#include "MathLib/MathTools.h"
12#include "NumLib/Exceptions.h"
13
14namespace MaterialLib
15{
16namespace Fracture
17{
18namespace Coulomb
19{
21{
22 double Kn = 0.0;
23 double Ks = 0.0;
24 double phi = 0.0; // friction angle
25 double psi = 0.0; // dilation angle
26 double c = 0.0;
27
28 template <typename MaterialProperties>
29 MaterialPropertyValues(MaterialProperties const& mp,
30 double const t,
32 {
33 Kn = mp.normal_stiffness(t, x)[0];
34 Ks = mp.shear_stiffness(t, x)[0];
35 auto constexpr degree =
36 boost::math::constants::degree<double>(); // pi/180
37 phi = mp.friction_angle(t, x)[0] * degree;
38 psi = mp.dilatancy_angle(t, x)[0] * degree;
39 c = mp.cohesion(t, x)[0];
40 }
41};
42
43template <int DisplacementDim>
45 double const t,
47 double const aperture0,
48 Eigen::Ref<Eigen::VectorXd const>
49 /*sigma0*/,
50 Eigen::Ref<Eigen::VectorXd const>
51 w_prev,
52 Eigen::Ref<Eigen::VectorXd const>
53 w,
54 Eigen::Ref<Eigen::VectorXd const>
55 sigma_prev,
56 Eigen::Ref<Eigen::VectorXd>
57 sigma,
58 Eigen::Ref<Eigen::MatrixXd>
59 Kep,
61 material_state_variables)
62{
63 assert(dynamic_cast<StateVariables<DisplacementDim> const*>(
64 &material_state_variables) != nullptr);
65
67 static_cast<StateVariables<DisplacementDim>&>(material_state_variables);
68
69 MaterialPropertyValues const mat(_mp, t, x);
70
71 const int index_ns = DisplacementDim - 1;
72 double const aperture = w[index_ns] + aperture0;
73
74 Eigen::MatrixXd Ke;
75 { // Elastic tangent stiffness
76 Ke = Eigen::MatrixXd::Zero(DisplacementDim, DisplacementDim);
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;
121 Kep(index_ns, index_ns) *= logPenaltyDerivative(
122 aperture0, aperture, _penalty_aperture_cutoff);
123 return;
124 }
125 }
126
127 auto yield_function_derivative = [&mat](Eigen::VectorXd const& s)
128 {
129 Eigen::Matrix<double, DisplacementDim, 1> dFs_dS;
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
137 auto plastic_potential_derivative = [&mat](Eigen::VectorXd const& s)
138 {
139 Eigen::Matrix<double, DisplacementDim, 1> dQs_dS;
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
148 Eigen::FullPivLU<Eigen::Matrix<double, 1, 1, Eigen::RowMajor>>
149 linear_solver;
150 using ResidualVectorType = Eigen::Matrix<double, 1, 1, Eigen::RowMajor>;
151 using JacobianMatrix = Eigen::Matrix<double, 1, 1, Eigen::RowMajor>;
152
153 JacobianMatrix jacobian;
154 ResidualVectorType solution;
155 solution << 0;
156
157 auto const update_residual = [&](ResidualVectorType& residual)
158 { residual[0] = yield_function(sigma); };
159
160 auto const update_jacobian = [&](JacobianMatrix& jacobian)
161 {
162 jacobian(0, 0) = -yield_function_derivative(sigma).transpose() *
163 Ke * plastic_potential_derivative(sigma);
164 };
165
166 auto const update_solution = [&](ResidualVectorType const& increment)
167 {
168 solution += increment;
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 +
173 solution[0] * plastic_potential_derivative(sigma);
174
175 sigma.noalias() = Ke * (w - w_prev - state.w_p + state.w_p_prev);
176
177 sigma.coeffRef(index_ns) *= logPenaltyDerivative(
178 aperture0, aperture, _penalty_aperture_cutoff);
179 sigma.noalias() += sigma_prev;
180 };
181
182 auto newton_solver = NumLib::NewtonRaphson(
183 linear_solver, update_jacobian, update_residual, update_solution,
185
186 auto const success_iterations = newton_solver.solve(jacobian);
187
188 if (!success_iterations)
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) *=
206 Eigen::RowVectorXd const A = yield_function_derivative(sigma).transpose() *
207 Ke /
208 (yield_function_derivative(sigma).transpose() *
209 Ke * plastic_potential_derivative(sigma));
210 Kep = Ke - Ke * plastic_potential_derivative(sigma) * A;
211}
212
213template class Coulomb<2>;
214template class Coulomb<3>;
215
216} // namespace Coulomb
217} // namespace Fracture
218} // namespace MaterialLib
NumLib::NewtonRaphsonSolverParameters const _nonlinear_solver_parameters
Definition Coulomb.h:129
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
Definition Coulomb.cpp:44
double logPenaltyDerivative(double const aperture0, double const aperture, double const aperture_cutoff)
Definition LogPenalty.h:17
MaterialPropertyValues(MaterialProperties const &mp, double const t, ParameterLib::SpatialPosition const &x)
Definition Coulomb.cpp:29
Eigen::Matrix< double, DisplacementDim, 1 > w_p
Definition Coulomb.h:29
Eigen::Matrix< double, DisplacementDim, 1 > w_p_prev
Definition Coulomb.h:33