OGS
MaterialPropertyLib::OrthotropicEmbeddedFracturePermeability< DisplacementDim > Class Template Referencefinal

Detailed Description

template<int DisplacementDim>
class MaterialPropertyLib::OrthotropicEmbeddedFracturePermeability< DisplacementDim >

Extended Permeability model based on Olivella&Alonso.

This property must be a medium property, it computes the permeability in dependence of the strain

The base model was proposed in [2] and it was further investigated in [30] . This extended Version features three orthotropic fracture planes.

The model takes the form of

\[ \mathbf{k} = k_\text{m} \mathbf{I} + \sum \limits_{i=1}^3 \frac{b_i}{a_i} \left( \frac{b_i^2}{12} - k_\text{m} \right) \left( \mathbf{I} - \mathbf{M}_i \right) \]

with

\[ \mathbf{M}_i = \vec{n}_i \otimes \vec{n}_i \]

and

\[ b_i = b_{i0} + \Delta b_i \\ \Delta b_i = a_i \langle \mathbf{\epsilon} : \mathbf{M}_i - \varepsilon_{0i} \rangle \]

where

\( k_\text{m} \) permeability of undisturbed material
\( b_i \) fracture aperture of each fracture plane
\( b_{i0} \) initial aperture of each fracture plane
\( a_i \) mean fracture distance of each fracture plane
\( \vec{n}_i \) fracture normal vector of each fracture plane
\( \varepsilon_{i0} \) threshold strain of each fracture plane

Definition at line 57 of file OrthotropicEmbeddedFracturePermeability.h.

#include <OrthotropicEmbeddedFracturePermeability.h>

Inheritance diagram for MaterialPropertyLib::OrthotropicEmbeddedFracturePermeability< DisplacementDim >:
[legend]
Collaboration diagram for MaterialPropertyLib::OrthotropicEmbeddedFracturePermeability< DisplacementDim >:
[legend]

Public Types

using SymmetricTensor
 

Public Member Functions

 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)
 
void checkScale () const override
 
PropertyDataType value (VariableArray const &variable_array, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const override
 
PropertyDataType dValue (VariableArray const &variable_array, Variable const variable, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const override
 
- Public Member Functions inherited from MaterialPropertyLib::Property
virtual ~Property ()
 
virtual PropertyDataType initialValue (ParameterLib::SpatialPosition const &pos, double const t) const
 
virtual PropertyDataType value () const
 
virtual PropertyDataType value (VariableArray const &variable_array, VariableArray const &variable_array_prev, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const
 
virtual PropertyDataType dValue (VariableArray const &variable_array, VariableArray const &variable_array_prev, Variable const variable, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const
 
virtual PropertyDataType d2Value (VariableArray const &variable_array, Variable const variable1, Variable const variable2, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const
 Default implementation: 2nd derivative of any constant property is zero.
 
virtual void setProperties (std::vector< std::unique_ptr< Phase > > const &phases)
 Default implementation:
 
void setScale (std::variant< Medium *, Phase *, Component * > scale)
 
template<typename T >
initialValue (ParameterLib::SpatialPosition const &pos, double const t) const
 
template<typename T >
value () const
 
template<typename T >
value (VariableArray const &variable_array, VariableArray const &variable_array_prev, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const
 
template<typename T >
value (VariableArray const &variable_array, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const
 
template<typename T >
dValue (VariableArray const &variable_array, VariableArray const &variable_array_prev, Variable const variable, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const
 
template<typename T >
dValue (VariableArray const &variable_array, Variable const variable, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const
 
template<typename T >
d2Value (VariableArray const &variable_array, Variable const &variable1, Variable const &variable2, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const
 

Private Attributes

Medium_medium = nullptr
 
std::vector< double > const _a
 
std::vector< double > const _e0
 
Eigen::Matrix< double, 3, 3 > const _n
 
ParameterLib::Parameter< double > const & _k
 
ParameterLib::Parameter< double > const & _phi_xy
 
ParameterLib::Parameter< double > const & _phi_yz
 
double const _jf
 

Additional Inherited Members

- Protected Attributes inherited from MaterialPropertyLib::Property
std::string name_
 
PropertyDataType value_
 The single value of a property.
 
PropertyDataType dvalue_
 
std::variant< Medium *, Phase *, Component * > scale_
 

Member Typedef Documentation

◆ SymmetricTensor

template<int DisplacementDim>
using MaterialPropertyLib::OrthotropicEmbeddedFracturePermeability< DisplacementDim >::SymmetricTensor
Initial value:
Eigen::Matrix<
double,
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.

Definition at line 80 of file OrthotropicEmbeddedFracturePermeability.h.

Constructor & Destructor Documentation

◆ OrthotropicEmbeddedFracturePermeability()

template<int DisplacementDim>
MaterialPropertyLib::OrthotropicEmbeddedFracturePermeability< DisplacementDim >::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 )

Definition at line 20 of file OrthotropicEmbeddedFracturePermeability.cpp.

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}

References MaterialPropertyLib::name, and MaterialPropertyLib::Property::name_.

Member Function Documentation

◆ checkScale()

template<int DisplacementDim>
void MaterialPropertyLib::OrthotropicEmbeddedFracturePermeability< DisplacementDim >::checkScale ( ) const
overridevirtual

Reimplemented from MaterialPropertyLib::Property.

Definition at line 42 of file OrthotropicEmbeddedFracturePermeability.cpp.

44{
45 if (!std::holds_alternative<Medium*>(scale_))
46 {
48 "The property 'OrthotropicEmbeddedFracturePermeability' is "
49 "implemented on the 'media' scale only.");
50 }
51}
#define OGS_FATAL(...)
Definition Error.h:26
std::variant< Medium *, Phase *, Component * > scale_
Definition Property.h:297

References OGS_FATAL.

◆ dValue()

template<int DisplacementDim>
PropertyDataType MaterialPropertyLib::OrthotropicEmbeddedFracturePermeability< DisplacementDim >::dValue ( VariableArray const & variable_array,
Variable const variable,
ParameterLib::SpatialPosition const & pos,
double const t,
double const dt ) const
overridevirtual

This virtual method will compute the property derivative value based on the variables that are passed as arguments with the default implementation using empty variables array for the previous time step.

The default implementation of this method only returns the property value derivative without altering it.

Reimplemented from MaterialPropertyLib::Property.

Definition at line 91 of file OrthotropicEmbeddedFracturePermeability.cpp.

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) *
126 MathLib::KelvinVector::tensorToKelvin<DisplacementDim>(
127 Eigen::Matrix3d::Identity() - M) *
128 MathLib::KelvinVector::tensorToKelvin<DisplacementDim>(M)
129 .transpose();
130 }
131 return Eigen::MatrixXd(_jf * result);
132}
PropertyDataType fromVector(std::vector< double > const &values)
Definition Property.cpp:23
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)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor > KelvinMatrixType

References MaterialPropertyLib::fromVector(), MathLib::KelvinVector::kelvinVectorToTensor(), MaterialPropertyLib::mechanical_strain, MaterialPropertyLib::VariableArray::mechanical_strain, and OGS_FATAL.

◆ value()

template<int DisplacementDim>
PropertyDataType MaterialPropertyLib::OrthotropicEmbeddedFracturePermeability< DisplacementDim >::value ( VariableArray const & variable_array,
ParameterLib::SpatialPosition const & pos,
double const t,
double const dt ) const
overridevirtual

This virtual method will compute the property value based on the variables that are passed as arguments with the default implementation using empty variables array for the previous time step.

Reimplemented from MaterialPropertyLib::Property.

Definition at line 55 of file OrthotropicEmbeddedFracturePermeability.cpp.

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}

References MaterialPropertyLib::fromVector(), MathLib::KelvinVector::kelvinVectorToTensor(), and MaterialPropertyLib::VariableArray::mechanical_strain.

Member Data Documentation

◆ _a

template<int DisplacementDim>
std::vector<double> const MaterialPropertyLib::OrthotropicEmbeddedFracturePermeability< DisplacementDim >::_a
private

Definition at line 61 of file OrthotropicEmbeddedFracturePermeability.h.

◆ _e0

template<int DisplacementDim>
std::vector<double> const MaterialPropertyLib::OrthotropicEmbeddedFracturePermeability< DisplacementDim >::_e0
private

Definition at line 62 of file OrthotropicEmbeddedFracturePermeability.h.

◆ _jf

template<int DisplacementDim>
double const MaterialPropertyLib::OrthotropicEmbeddedFracturePermeability< DisplacementDim >::_jf
private

Definition at line 67 of file OrthotropicEmbeddedFracturePermeability.h.

◆ _k

template<int DisplacementDim>
ParameterLib::Parameter<double> const& MaterialPropertyLib::OrthotropicEmbeddedFracturePermeability< DisplacementDim >::_k
private

Definition at line 64 of file OrthotropicEmbeddedFracturePermeability.h.

◆ _medium

template<int DisplacementDim>
Medium* MaterialPropertyLib::OrthotropicEmbeddedFracturePermeability< DisplacementDim >::_medium = nullptr
private

Definition at line 60 of file OrthotropicEmbeddedFracturePermeability.h.

◆ _n

template<int DisplacementDim>
Eigen::Matrix<double, 3, 3> const MaterialPropertyLib::OrthotropicEmbeddedFracturePermeability< DisplacementDim >::_n
private

Definition at line 63 of file OrthotropicEmbeddedFracturePermeability.h.

◆ _phi_xy

template<int DisplacementDim>
ParameterLib::Parameter<double> const& MaterialPropertyLib::OrthotropicEmbeddedFracturePermeability< DisplacementDim >::_phi_xy
private

Definition at line 65 of file OrthotropicEmbeddedFracturePermeability.h.

◆ _phi_yz

template<int DisplacementDim>
ParameterLib::Parameter<double> const& MaterialPropertyLib::OrthotropicEmbeddedFracturePermeability< DisplacementDim >::_phi_yz
private

Definition at line 66 of file OrthotropicEmbeddedFracturePermeability.h.


The documentation for this class was generated from the following files: