Loading [MathJax]/extensions/tex2jax.js
OGS
OrthotropicEmbeddedFracturePermeability.cpp
Go to the documentation of this file.
1 
12 
13 #include "MaterialLib/MPL/Medium.h"
15 
16 namespace MaterialPropertyLib
17 {
18 template <int DisplacementDim>
21  std::string name,
22  std::vector<double> const& mean_fracture_distances,
23  std::vector<double> const& threshold_strains,
24  Eigen::Matrix<double, 3, 3> const fracture_normals,
25  ParameterLib::Parameter<double> const& intrinsic_permeability,
26  ParameterLib::Parameter<double> const& fracture_rotation_xy,
27  ParameterLib::Parameter<double> const& fracture_rotation_yz)
28  : _a(mean_fracture_distances),
29  _e0(threshold_strains),
30  _n(fracture_normals),
31  _k(intrinsic_permeability),
32  _phi_xy(fracture_rotation_xy),
33  _phi_yz(fracture_rotation_yz)
34 {
35  name_ = std::move(name);
36 }
37 
38 template <int DisplacementDim>
40  const
41 {
42  if (!std::holds_alternative<Medium*>(scale_))
43  {
44  OGS_FATAL(
45  "The property 'OrthotropicEmbeddedFracturePermeability' is "
46  "implemented on the 'media' scale only.");
47  }
48 }
49 
50 template <int DisplacementDim>
53  VariableArray const& variable_array,
54  ParameterLib::SpatialPosition const& pos, double const t,
55  double const /*dt*/) const
56 {
59  variable_array[static_cast<int>(Variable::mechanical_strain)]));
60  double const k = std::get<double>(fromVector(_k(t, pos)));
61  double const _b0 = std::sqrt(12.0 * k);
62 
63  double const phi_xy = std::get<double>(fromVector(_phi_xy(t, pos)));
64  double const phi_yz = std::get<double>(fromVector(_phi_yz(t, pos)));
65  auto const rotMat_xy = Eigen::AngleAxisd(phi_xy, Eigen::Vector3d::UnitZ());
66  auto const rotMat_yz = Eigen::AngleAxisd(phi_yz, Eigen::Vector3d::UnitX());
67 
68  Eigen::Matrix3d result = Eigen::Vector3d::Constant(k).asDiagonal();
69 
70  for (int i = 0; i < 3; i++)
71  {
72  Eigen::Vector3d const n_i = rotMat_yz * (rotMat_xy * _n.col(i));
73  double const e_n = (eps * n_i).dot(n_i.transpose());
74  double const H_de = (e_n > _e0[i]) ? 1.0 : 0.0;
75  double const b_f = _b0 + H_de * _a[i] * (e_n - _e0[i]);
76 
77  result += H_de * (b_f / _a[i]) * ((b_f * b_f / 12.0) - k) *
78  (Eigen::Matrix3d::Identity() - n_i * n_i.transpose());
79  }
80 
81  return result.template topLeftCorner<DisplacementDim, DisplacementDim>()
82  .eval();
83 }
84 template <int DisplacementDim>
87  VariableArray const& variable_array, Variable const primary_variable,
88  ParameterLib::SpatialPosition const& pos, double const t,
89  double const /*dt*/) const
90 {
91  if (primary_variable != Variable::mechanical_strain)
92  {
93  OGS_FATAL(
94  "OrthotropicEmbeddedFracturePermeability::dValue is implemented "
95  "for derivatives with respect to strain only.");
96  }
97 
98  auto const eps = formEigenTensor<3>(std::get<SymmetricTensor>(
99  variable_array[static_cast<int>(Variable::mechanical_strain)]));
100  double const k = std::get<double>(fromVector(_k(t, pos)));
101  double const _b0 = std::sqrt(12.0 * k);
102 
103  double const phi_xy = std::get<double>(fromVector(_phi_xy(t, pos)));
104  double const phi_yz = std::get<double>(fromVector(_phi_yz(t, pos)));
105  auto const rotMat_xy = Eigen::AngleAxisd(phi_xy, Eigen::Vector3d::UnitZ());
106  auto const rotMat_yz = Eigen::AngleAxisd(phi_yz, Eigen::Vector3d::UnitX());
107 
108  Eigen::Matrix3d result = Eigen::Matrix3d::Zero();
109 
110  for (int i = 0; i < 3; i++)
111  {
112  Eigen::Vector3d const n_i = rotMat_yz * (rotMat_xy * _n.col(i));
113  Eigen::Matrix3d const M = n_i * n_i.transpose();
114  double const e_n = (eps * n_i).dot(n_i.transpose());
115  double const H_de = (e_n > _e0[i]) ? 1.0 : 0.0;
116  double const b_f = _b0 + H_de * _a[i] * (e_n - _e0[i]);
117 
118  result += H_de * (b_f * b_f / 4.0 - k) *
119  (Eigen::Matrix3d::Identity() - M) * M;
120  }
121  return result.template topLeftCorner<DisplacementDim, DisplacementDim>()
122  .eval();
123 }
124 
127 } // namespace MaterialPropertyLib
#define OGS_FATAL(...)
Definition: Error.h:26
PropertyDataType dValue(VariableArray const &variable_array, Variable const primary_variable, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const override
OrthotropicEmbeddedFracturePermeability(std::string name, std::vector< double > const &mean_fracture_distances, std::vector< double > const &threshold_strains, Eigen::Matrix< double, 3, 3 > const fracture_normals, ParameterLib::Parameter< double > const &intrinsic_permeability, ParameterLib::Parameter< double > const &fracture_rotation_xy, ParameterLib::Parameter< double > const &fracture_rotation_yz)
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
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)