OGS
|
Permeability model as proposed by Olivella&Alonso.
This property must be a medium property, it computes the permeability in dependence of the strain
The model was proposed in [2] and it was further investigated in [31] .
The model takes the form of
\[ \mathbf{k} = k_\text{m} \mathbf{I} + \frac{b}{a} \left( \frac{b^2}{12} - k_\text{m} \right) \left( \mathbf{I} - \mathbf{M} \right) \]
with
\[ \mathbf{M} = \vec{n} \otimes \vec{n} \]
and
\[ b = b_0 + \Delta b \\ \Delta b = a \langle \varepsilon_\text{n}- \varepsilon_0 \rangle \]
where
\( k_\text{m} \) | permeability of undisturbed material |
\( b \) | fracture aperture |
\( b_0 \) | initial aperture |
\( a \) | mean fracture distance |
\( \vec{n} \) | fracture normal vector |
\( \varepsilon_n \) | strain in fracture normal direction |
\( \varepsilon_0 \) | threshold strain |
Definition at line 50 of file EmbeddedFracturePermeability.h.
#include <EmbeddedFracturePermeability.h>
Public Types | |
using | SymmetricTensor = Eigen::Matrix<double, KelvinVectorSize, 1> |
Public Member Functions | |
EmbeddedFracturePermeability (std::string name, Eigen::Matrix< double, 3, 1 > const fracture_normal, bool const fracture_normal_is_constant, double const intrinsic_permeability, double const initial_aperture, double const mean_fracture_distance, double const threshold_strain, 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 |
Static Public Attributes | |
static int const | KelvinVectorSize |
Private Attributes | |
Medium * | _medium = nullptr |
Eigen::Matrix< double, 3, 1 > const | _n |
bool const | _n_const |
double const | _k |
double const | _b0 |
double const | _a |
double const | _e0 |
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::EmbeddedFracturePermeability< DisplacementDim >::SymmetricTensor = Eigen::Matrix<double, KelvinVectorSize, 1> |
Definition at line 79 of file EmbeddedFracturePermeability.h.
MaterialPropertyLib::EmbeddedFracturePermeability< DisplacementDim >::EmbeddedFracturePermeability | ( | std::string | name, |
Eigen::Matrix< double, 3, 1 > const | fracture_normal, | ||
bool const | fracture_normal_is_constant, | ||
double const | intrinsic_permeability, | ||
double const | initial_aperture, | ||
double const | mean_fracture_distance, | ||
double const | threshold_strain, | ||
ParameterLib::Parameter< double > const & | fracture_rotation_xy, | ||
ParameterLib::Parameter< double > const & | fracture_rotation_yz, | ||
double const | jacobian_factor ) |
Definition at line 22 of file EmbeddedFracturePermeability.cpp.
References MaterialPropertyLib::name, and MaterialPropertyLib::Property::name_.
|
overridevirtual |
Reimplemented from MaterialPropertyLib::Property.
Definition at line 47 of file EmbeddedFracturePermeability.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 97 of file EmbeddedFracturePermeability.cpp.
References MaterialPropertyLib::formEigenTensor< 3 >(), MathLib::KelvinVector::kelvinVectorToTensor(), MaterialPropertyLib::mechanical_strain, MaterialPropertyLib::VariableArray::mechanical_strain, OGS_FATAL, MathLib::KelvinVector::tensorToKelvin(), and MaterialPropertyLib::VariableArray::total_stress.
|
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 58 of file EmbeddedFracturePermeability.cpp.
References MaterialPropertyLib::formEigenTensor< 3 >(), MathLib::KelvinVector::kelvinVectorToTensor(), MaterialPropertyLib::VariableArray::mechanical_strain, and MaterialPropertyLib::VariableArray::total_stress.
|
private |
Definition at line 58 of file EmbeddedFracturePermeability.h.
|
private |
Definition at line 57 of file EmbeddedFracturePermeability.h.
|
private |
Definition at line 59 of file EmbeddedFracturePermeability.h.
|
private |
Definition at line 62 of file EmbeddedFracturePermeability.h.
|
private |
Definition at line 56 of file EmbeddedFracturePermeability.h.
|
private |
Definition at line 53 of file EmbeddedFracturePermeability.h.
|
private |
Definition at line 54 of file EmbeddedFracturePermeability.h.
|
private |
Definition at line 55 of file EmbeddedFracturePermeability.h.
|
private |
Definition at line 60 of file EmbeddedFracturePermeability.h.
|
private |
Definition at line 61 of file EmbeddedFracturePermeability.h.
|
static |
Definition at line 77 of file EmbeddedFracturePermeability.h.