OGS
PermeabilityMohrCoulombFailureIndexModel.cpp
Go to the documentation of this file.
1 
13 
14 #include <boost/math/constants/constants.hpp>
15 #include <cmath>
16 #include <limits>
17 
18 #include "BaseLib/Error.h"
19 #include "MaterialLib/MPL/Medium.h"
22 #include "MathLib/KelvinVector.h"
23 #include "MathLib/MathTools.h"
25 #include "ParameterLib/Parameter.h"
26 
27 namespace MaterialPropertyLib
28 {
29 template <int DisplacementDim>
32  std::string name, ParameterLib::Parameter<double> const& k0,
33  double const kr, double const b, double const c, double const phi,
34  double const k_max, double const t_sigma_max,
35  ParameterLib::CoordinateSystem const* const local_coordinate_system)
36  : k0_(k0),
37  kr_(kr),
38  b_(b),
39  c_(c),
40  phi_(boost::math::constants::degree<double>() * phi),
41  k_max_(k_max),
42  t_sigma_max_(t_sigma_max),
43  local_coordinate_system_(local_coordinate_system)
44 {
45  const double t_sigma_upper = c_ / std::tan(phi_);
46  if (t_sigma_max_ <= 0.0 || t_sigma_max_ > t_sigma_upper ||
47  std::fabs(t_sigma_max_ - t_sigma_upper) <
48  std::numeric_limits<double>::epsilon())
49  {
50  OGS_FATAL(
51  "Tensile strength parameter of {:e} is out of the range (0, "
52  "c/tan(phi)) = (0, {:e})",
53  t_sigma_max_, t_sigma_upper);
54  }
55 
56  name_ = std::move(name);
57 }
58 
59 template <int DisplacementDim>
61  const
62 {
63  if (!std::holds_alternative<Medium*>(scale_))
64  {
65  OGS_FATAL(
66  "The property 'PermeabilityMohrCoulombFailureIndexModel' is "
67  "implemented on the 'medium' scale only.");
68  }
69 }
70 
71 template <int DisplacementDim>
74  VariableArray const& variable_array,
75  ParameterLib::SpatialPosition const& pos, double const t,
76  double const /*dt*/) const
77 {
78  auto const& stress_vector = std::get<SymmetricTensor<DisplacementDim>>(
79  variable_array[static_cast<int>(Variable::total_stress)]);
80 
81  auto const& stress_tensor =
82  formEigenTensor<3>(static_cast<PropertyDataType>(stress_vector));
83 
84  Eigen::SelfAdjointEigenSolver<Eigen::Matrix<double, 3, 3>>
85  eigenvalue_solver(stress_tensor);
86 
87  // Principle stress
88  auto const sigma = eigenvalue_solver.eigenvalues();
89 
90  auto k_data = k0_(t, pos);
91 
92  double const max_sigma = std::max(std::fabs(sigma[0]), std::fabs(sigma[2]));
93 
94  if (max_sigma < std::numeric_limits<double>::epsilon())
95  {
96  return fromVector(k_data);
97  }
98 
99  double const sigma_m = 0.5 * (sigma[2] + sigma[0]);
100 
101  double const tau_m = 0.5 * std::fabs(sigma[2] - sigma[0]);
102  double f = 0.0;
103  if (sigma_m > t_sigma_max_)
104  {
105  // tensile failure criterion
106  f = sigma_m / t_sigma_max_;
107 
108  double const tau_tt =
109  c_ * std::cos(phi_) - t_sigma_max_ * std::sin(phi_);
110 
111  f = std::max(f, tau_m / tau_tt);
112  }
113  else
114  {
115  // Mohr Coulomb failure criterion
116  f = tau_m / (c_ * std::cos(phi_) - sigma_m * std::sin(phi_));
117  }
118 
119  if (f >= 1.0)
120  {
121  const double exp_value = std::exp(b_ * f);
122  for (auto& k_i : k_data)
123  {
124  k_i = std::min(k_i + kr_ * exp_value, k_max_);
125  }
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 
150 template <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
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:72
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 > > PropertyDataType
Definition: Property.h:35
template Eigen::Matrix< double, 3, 3 > formEigenTensor< 3 >(MaterialPropertyLib::PropertyDataType const &values)
std::array< VariableType, static_cast< int >(Variable::number_of_variables)> VariableArray
Definition: VariableType.h:108