OGS
|
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 [31] . 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>
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 > | |
T | initialValue (ParameterLib::SpatialPosition const &pos, double const t) const |
template<typename T > | |
T | value () const |
template<typename T > | |
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 > | |
T | value (VariableArray const &variable_array, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const |
template<typename T > | |
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 > | |
T | dValue (VariableArray const &variable_array, Variable const variable, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const |
template<typename T > | |
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_ |
using MaterialPropertyLib::OrthotropicEmbeddedFracturePermeability< DisplacementDim >::SymmetricTensor |
Definition at line 80 of file OrthotropicEmbeddedFracturePermeability.h.
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.
References MaterialPropertyLib::name, and MaterialPropertyLib::Property::name_.
|
overridevirtual |
Reimplemented from MaterialPropertyLib::Property.
Definition at line 42 of file OrthotropicEmbeddedFracturePermeability.cpp.
References OGS_FATAL.
|
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.
References MaterialPropertyLib::fromVector(), MathLib::KelvinVector::kelvinVectorToTensor(), MaterialPropertyLib::mechanical_strain, MaterialPropertyLib::VariableArray::mechanical_strain, and OGS_FATAL.
|
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.
References MaterialPropertyLib::fromVector(), MathLib::KelvinVector::kelvinVectorToTensor(), and MaterialPropertyLib::VariableArray::mechanical_strain.
|
private |
Definition at line 61 of file OrthotropicEmbeddedFracturePermeability.h.
|
private |
Definition at line 62 of file OrthotropicEmbeddedFracturePermeability.h.
|
private |
Definition at line 67 of file OrthotropicEmbeddedFracturePermeability.h.
|
private |
Definition at line 64 of file OrthotropicEmbeddedFracturePermeability.h.
|
private |
Definition at line 60 of file OrthotropicEmbeddedFracturePermeability.h.
|
private |
Definition at line 63 of file OrthotropicEmbeddedFracturePermeability.h.
|
private |
Definition at line 65 of file OrthotropicEmbeddedFracturePermeability.h.
|
private |
Definition at line 66 of file OrthotropicEmbeddedFracturePermeability.h.