OGS
EmbeddedFracturePermeability.cpp
Go to the documentation of this file.
1 
12 
13 #include "MaterialLib/MPL/Medium.h"
15 #include "MathLib/KelvinVector.h"
16 
17 namespace MaterialPropertyLib
18 {
19 template <int DisplacementDim>
21  std::string name,
22  Eigen::Matrix<double, 3, 1> const fracture_normal,
23  bool const fracture_normal_is_constant,
24  double const intrinsic_permeability,
25  double const initial_aperture,
26  double const mean_fracture_distance,
27  double const threshold_strain)
28  : _n(fracture_normal),
29  _n_const(fracture_normal_is_constant),
30  _k(intrinsic_permeability),
31  _b0(initial_aperture),
32  _a(mean_fracture_distance),
33  _e0(threshold_strain)
34 {
35  name_ = std::move(name);
36 }
37 
38 template <int DisplacementDim>
40 {
41  if (!std::holds_alternative<Medium*>(scale_))
42  {
43  OGS_FATAL(
44  "The property 'EmbeddedFracturePermeability' is "
45  "implemented on the 'medium' scale only.");
46  }
47 }
48 
49 template <int DisplacementDim>
51  VariableArray const& variable_array,
52  ParameterLib::SpatialPosition const& /*pos*/, double const /*t*/,
53  double const /*dt*/) const
54 {
55  Eigen::Matrix<double, 3, 1> const n = [&]
56  {
57  if (_n_const)
58  {
59  return _n;
60  }
61  auto const sigma = formEigenTensor<3>(std::get<SymmetricTensor>(
62  variable_array[static_cast<int>(Variable::total_stress)]));
63  Eigen::SelfAdjointEigenSolver<Eigen::Matrix<double, 3, 3>> e_s(sigma);
64  return (Eigen::Matrix<double, 3, 1>)e_s.eigenvectors().col(2);
65  }();
66 
69  variable_array[static_cast<int>(Variable::mechanical_strain)]));
70  double const e_n = (eps * n).dot(n.transpose());
71  double const H_de = (e_n > _e0) ? 1.0 : 0.0;
72  double const b_f = _b0 + H_de * _a * (e_n - _e0);
73  double const coeff = H_de * (b_f / _a) * ((b_f * b_f / 12.0) - _k);
74 
75  Eigen::Matrix3d I = Eigen::Matrix3d::Identity();
76  return (coeff * (I - n * n.transpose()) + _k * I)
77  .template topLeftCorner<DisplacementDim, DisplacementDim>()
78  .eval();
79 }
80 
81 template <int DisplacementDim>
83  VariableArray const& variable_array, Variable const primary_variable,
84  ParameterLib::SpatialPosition const& /*pos*/, double const /*t*/,
85  double const /*dt*/) const
86 {
87  if (primary_variable != Variable::mechanical_strain)
88  {
89  OGS_FATAL(
90  "EmbeddedFracturePermeability::dValue is implemented for "
91  "derivatives with respect to strain only.");
92  }
93 
94  Eigen::Matrix<double, 3, 1> const n = [&]
95  {
96  if (_n_const)
97  {
98  return _n;
99  }
100  auto const sigma = formEigenTensor<3>(std::get<SymmetricTensor>(
101  variable_array[static_cast<int>(Variable::total_stress)]));
102  Eigen::SelfAdjointEigenSolver<Eigen::Matrix<double, 3, 3>> e_s(sigma);
103  return (Eigen::Matrix<double, 3, 1>)e_s.eigenvectors().col(2);
104  }();
105 
108  variable_array[static_cast<int>(Variable::mechanical_strain)]));
109  double const e_n = (eps * n).dot(n.transpose());
110  double const H_de = (e_n > _e0) ? 1.0 : 0.0;
111  double const b_f = _b0 + H_de * _a * (e_n - _e0);
112 
113  Eigen::Matrix3d const M = n * n.transpose();
114  return (H_de * (b_f * b_f / 4 - _k) * (Eigen::Matrix3d::Identity() - M) * M)
115  .template topLeftCorner<DisplacementDim, DisplacementDim>()
116  .eval();
117 }
118 
119 template class EmbeddedFracturePermeability<2>;
120 template class EmbeddedFracturePermeability<3>;
121 } // namespace MaterialPropertyLib
#define OGS_FATAL(...)
Definition: Error.h:26
Permeability model as proposed by Olivella&Alonso.
PropertyDataType dValue(VariableArray const &variable_array, Variable const primary_variable, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const override
EmbeddedFracturePermeability(std::string name, Eigen::Matrix< double, 3, 1 > const fracture_normal, bool const fracture_normal_is_constant, double const intrinsic_permeability, double const initial_aperture, double const mean_fracture_distance, double const threshold_strain)
virtual PropertyDataType value() const
Definition: Property.cpp:72
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
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Definition: KelvinVector.h:48
Eigen::Matrix< double, 3, 3 > kelvinVectorToTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)