OGS
PermeabilityMohrCoulombFailureIndexModel.cpp
Go to the documentation of this file.
1
13
14#include <Eigen/Eigenvalues>
15#include <boost/math/constants/constants.hpp>
16#include <cmath>
17#include <limits>
18
19#include "BaseLib/Error.h"
24#include "MathLib/MathTools.h"
27
28namespace MaterialPropertyLib
29{
30template <int DisplacementDim>
33 std::string name, ParameterLib::Parameter<double> const& k0,
34 double const kr, double const b, double const c, double const phi,
35 double const k_max, double const t_sigma_max,
36 ParameterLib::CoordinateSystem const* const local_coordinate_system)
37 : k0_(k0),
38 kr_(kr),
39 b_(b),
40 c_(c),
41 phi_(boost::math::constants::degree<double>() * phi),
42 k_max_(k_max),
43 t_sigma_max_(t_sigma_max),
44 local_coordinate_system_(local_coordinate_system)
45{
46 const double t_sigma_upper = c_ / std::tan(phi_);
47 if (t_sigma_max_ <= 0.0 || t_sigma_max_ > t_sigma_upper ||
48 std::fabs(t_sigma_max_ - t_sigma_upper) <
49 std::numeric_limits<double>::epsilon())
50 {
52 "Tensile strength parameter of {:e} is out of the range (0, "
53 "c/tan(phi)) = (0, {:e})",
54 t_sigma_max_, t_sigma_upper);
55 }
56
57 name_ = std::move(name);
58}
59
60template <int DisplacementDim>
62 const
63{
64 if (!std::holds_alternative<Medium*>(scale_))
65 {
67 "The property 'PermeabilityMohrCoulombFailureIndexModel' is "
68 "implemented on the 'medium' scale only.");
69 }
70}
71
72template <int DisplacementDim>
75 VariableArray const& variable_array,
76 ParameterLib::SpatialPosition const& pos, double const t,
77 double const /*dt*/) const
78{
79 auto const& stress_vector =
80 std::get<SymmetricTensor<DisplacementDim>>(variable_array.total_stress);
81
82 auto const& stress_tensor =
83 formEigenTensor<3>(static_cast<PropertyDataType>(stress_vector));
84
85 Eigen::SelfAdjointEigenSolver<Eigen::Matrix<double, 3, 3>>
86 eigenvalue_solver(stress_tensor);
87
88 // Principle stress
89 auto const sigma = eigenvalue_solver.eigenvalues();
90
91 auto k_data = k0_(t, pos);
92
93 double const max_sigma = std::max(std::fabs(sigma[0]), std::fabs(sigma[2]));
94
95 if (max_sigma < std::numeric_limits<double>::epsilon())
96 {
97 return fromVector(k_data);
98 }
99
100 double const sigma_m = 0.5 * (sigma[2] + sigma[0]);
101
102 double const tau_m = 0.5 * std::fabs(sigma[2] - sigma[0]);
103 double f = 0.0;
104 if (sigma_m > t_sigma_max_)
105 {
106 // tensile failure criterion
107 f = sigma_m / t_sigma_max_;
108
109 double const tau_tt =
110 c_ * std::cos(phi_) - t_sigma_max_ * std::sin(phi_);
111
112 f = std::max(f, tau_m / tau_tt);
113 }
114 else
115 {
116 // Mohr Coulomb failure criterion
117 f = tau_m / (c_ * std::cos(phi_) - sigma_m * std::sin(phi_));
118 }
119
120 if (f >= 1.0)
121 {
122 const double exp_value = std::exp(b_ * f);
123 std::transform(begin(k_data), end(k_data), begin(k_data),
124 [&](double const k_i)
125 { return std::min(k_i + kr_ * exp_value, k_max_); });
126 }
127
128 // Local coordinate transformation is only applied for the case that the
129 // initial intrinsic permeability is given with orthotropic assumption.
130 if (local_coordinate_system_ && (k_data.size() == DisplacementDim))
131 {
132 Eigen::Matrix<double, DisplacementDim, DisplacementDim> const e =
133 local_coordinate_system_->transformation<DisplacementDim>(pos);
134 Eigen::Matrix<double, DisplacementDim, DisplacementDim> k =
135 Eigen::Matrix<double, DisplacementDim, DisplacementDim>::Zero();
136
137 for (int i = 0; i < DisplacementDim; ++i)
138 {
139 Eigen::Matrix<double, DisplacementDim, DisplacementDim> const
140 ei_otimes_ei = e.col(i) * e.col(i).transpose();
141
142 k += k_data[i] * ei_otimes_ei;
143 }
144 return k;
145 }
146
147 return fromVector(k_data);
148}
149
150template <int DisplacementDim>
153 VariableArray const& /*variable_array*/, Variable const variable,
154 ParameterLib::SpatialPosition const& /*pos*/, double const /*t*/,
155 double const /*dt*/) const
156{
157 if (variable == Variable::mechanical_strain)
158 {
159 return 0.;
160 }
161
162 OGS_FATAL(
163 "The derivative of the intrinsic permeability k(sigma, ...) with "
164 "respect to stress tensor (sigma) is not implemented because that "
165 "dk/du is normally omitted.");
166}
169} // namespace MaterialPropertyLib
#define OGS_FATAL(...)
Definition Error.h:26
A failure index dependent permeability model wangetalperm2020.
PermeabilityMohrCoulombFailureIndexModel(std::string name, ParameterLib::Parameter< double > const &k0, double const kr, double const b, double const c, double const phi, double const k_max, double const t_sigma_max, ParameterLib::CoordinateSystem const *const local_coordinate_system)
PropertyDataType dValue(VariableArray const &variable_array, Variable const variable, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const override
virtual PropertyDataType value() const
Definition Property.cpp:76
template Eigen::Matrix< double, 3, 3 > formEigenTensor< 3 >(MaterialPropertyLib::PropertyDataType const &values)
PropertyDataType fromVector(std::vector< double > const &values)
Definition Property.cpp:23
std::variant< double, Eigen::Matrix< double, 2, 1 >, Eigen::Matrix< double, 3, 1 >, Eigen::Matrix< double, 2, 2 >, Eigen::Matrix< double, 3, 3 >, Eigen::Matrix< double, 4, 1 >, Eigen::Matrix< double, 6, 1 >, Eigen::MatrixXd > PropertyDataType
Definition Property.h:31
A local coordinate system used for tensor transformations.