18 template <
int DisplacementDim>
22 std::vector<double>
const& mean_fracture_distances,
23 std::vector<double>
const& threshold_strains,
24 Eigen::Matrix<double, 3, 3>
const fracture_normals,
28 : _a(mean_fracture_distances),
29 _e0(threshold_strains),
31 _k(intrinsic_permeability),
32 _phi_xy(fracture_rotation_xy),
33 _phi_yz(fracture_rotation_yz)
38 template <
int DisplacementDim>
42 if (!std::holds_alternative<Medium*>(scale_))
45 "The property 'OrthotropicEmbeddedFracturePermeability' is "
46 "implemented on the 'media' scale only.");
50 template <
int DisplacementDim>
60 double const k = std::get<double>(
fromVector(_k(t, pos)));
61 double const _b0 = std::sqrt(12.0 * k);
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());
68 Eigen::Matrix3d result = Eigen::Vector3d::Constant(k).asDiagonal();
70 for (
int i = 0; i < 3; i++)
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]);
77 result += H_de * (b_f / _a[i]) * ((b_f * b_f / 12.0) - k) *
78 (Eigen::Matrix3d::Identity() - n_i * n_i.transpose());
81 return result.template topLeftCorner<DisplacementDim, DisplacementDim>()
84 template <
int DisplacementDim>
94 "OrthotropicEmbeddedFracturePermeability::dValue is implemented "
95 "for derivatives with respect to strain only.");
100 double const k = std::get<double>(
fromVector(_k(t, pos)));
101 double const _b0 = std::sqrt(12.0 * k);
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());
108 Eigen::Matrix3d result = Eigen::Matrix3d::Zero();
110 for (
int i = 0; i < 3; i++)
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]);
118 result += H_de * (b_f * b_f / 4.0 - k) *
119 (Eigen::Matrix3d::Identity() - M) * M;
121 return result.template topLeftCorner<DisplacementDim, DisplacementDim>()
Extended Permeability model based on Olivella&Alonso.
PropertyDataType dValue(VariableArray const &variable_array, Variable const primary_variable, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const override
void checkScale() 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
PropertyDataType fromVector(std::vector< double > const &values)
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
template Eigen::Matrix< double, 3, 3 > formEigenTensor< 3 >(MaterialPropertyLib::PropertyDataType const &values)
std::array< VariableType, static_cast< int >(Variable::number_of_variables)> VariableArray
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)