OGS
OrthotropicEmbeddedFracturePermeability.cpp
Go to the documentation of this file.
1
12
13#include <Eigen/Geometry>
14
16
17namespace MaterialPropertyLib
18{
19template <int DisplacementDim>
22 std::string name,
23 std::vector<double> const& mean_fracture_distances,
24 std::vector<double> const& threshold_strains,
25 Eigen::Matrix<double, 3, 3> const fracture_normals,
26 ParameterLib::Parameter<double> const& intrinsic_permeability,
27 ParameterLib::Parameter<double> const& fracture_rotation_xy,
28 ParameterLib::Parameter<double> const& fracture_rotation_yz,
29 double const jacobian_factor)
30 : _a(mean_fracture_distances),
31 _e0(threshold_strains),
32 _n(fracture_normals),
33 _k(intrinsic_permeability),
34 _phi_xy(fracture_rotation_xy),
35 _phi_yz(fracture_rotation_yz),
36 _jf(jacobian_factor)
37{
38 name_ = std::move(name);
39}
40
41template <int DisplacementDim>
43 const
44{
45 if (!std::holds_alternative<Medium*>(scale_))
46 {
48 "The property 'OrthotropicEmbeddedFracturePermeability' is "
49 "implemented on the 'media' scale only.");
50 }
51}
52
53template <int DisplacementDim>
56 VariableArray const& variable_array,
57 ParameterLib::SpatialPosition const& pos, double const t,
58 double const /*dt*/) const
59{
62 variable_array.mechanical_strain));
63 double const k = std::get<double>(fromVector(_k(t, pos)));
64 double const _b0 = std::sqrt(12.0 * k);
65
66 double const phi_xy = std::get<double>(fromVector(_phi_xy(t, pos)));
67 double const phi_yz = std::get<double>(fromVector(_phi_yz(t, pos)));
68 auto const rotMat_xy = Eigen::AngleAxisd(phi_xy, Eigen::Vector3d::UnitZ());
69 auto const rotMat_yz = Eigen::AngleAxisd(phi_yz, Eigen::Vector3d::UnitX());
70
71 Eigen::Matrix3d result = Eigen::Vector3d::Constant(k).asDiagonal();
72
73 for (int i = 0; i < 3; i++)
74 {
75 Eigen::Vector3d const n_i = rotMat_yz * (rotMat_xy * _n.col(i));
76 double const e_n = (eps * n_i).dot(n_i.transpose());
77 double const H_de = (e_n > _e0[i]) ? 1.0 : 0.0;
78 double const b_f = _b0 + H_de * _a[i] * (e_n - _e0[i]);
79
80 // The H_de factor is only valid as long as _b0 equals
81 // std::sqrt(12.0 * k), else it has to be dropped.
82 result += H_de * (b_f / _a[i]) * ((b_f * b_f / 12.0) - k) *
83 (Eigen::Matrix3d::Identity() - n_i * n_i.transpose());
84 }
85
86 return result.template topLeftCorner<DisplacementDim, DisplacementDim>()
87 .eval();
88}
89template <int DisplacementDim>
92 VariableArray const& variable_array, Variable const variable,
93 ParameterLib::SpatialPosition const& pos, double const t,
94 double const /*dt*/) const
95{
96 if (variable != Variable::mechanical_strain)
97 {
99 "OrthotropicEmbeddedFracturePermeability::dValue is implemented "
100 "for derivatives with respect to strain only.");
101 }
102
105 variable_array.mechanical_strain));
106 double const k = std::get<double>(fromVector(_k(t, pos)));
107 double const _b0 = std::sqrt(12.0 * k);
108
109 double const phi_xy = std::get<double>(fromVector(_phi_xy(t, pos)));
110 double const phi_yz = std::get<double>(fromVector(_phi_yz(t, pos)));
111 auto const rotMat_xy = Eigen::AngleAxisd(phi_xy, Eigen::Vector3d::UnitZ());
112 auto const rotMat_yz = Eigen::AngleAxisd(phi_yz, Eigen::Vector3d::UnitX());
113
116
117 for (int i = 0; i < 3; i++)
118 {
119 Eigen::Vector3d const n_i = rotMat_yz * (rotMat_xy * _n.col(i));
120 Eigen::Matrix3d const M = n_i * n_i.transpose();
121 double const e_n = (eps * n_i).dot(n_i.transpose());
122 double const H_de = (e_n > _e0[i]) ? 1.0 : 0.0;
123 double const b_f = _b0 + H_de * _a[i] * (e_n - _e0[i]);
124
125 result += H_de * (b_f * b_f / 4.0 - k) *
127 Eigen::Matrix3d::Identity() - M) *
129 .transpose();
130 }
131 return Eigen::MatrixXd(_jf * result);
132}
133
136} // namespace MaterialPropertyLib
#define OGS_FATAL(...)
Definition Error.h:26
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, double const jacobian_factor)
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
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
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Eigen::Matrix< double, 3, 3 > kelvinVectorToTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
KelvinVectorType< DisplacementDim > tensorToKelvin(Eigen::Matrix< double, 3, 3 > const &m)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor > KelvinMatrixType