OGS
|
A failure index dependent permeability model [42].
\[ \mathbf{k} = \mathbf{k}_0+ H(f-1) k_\text{r} \mathrm{e}^{b f}\mathbf{I}\]
where \(\mathbf{k}_0\) is the intrinsic permeability of the undamaged material, \(H\) is the Heaviside step function, \(f\) is the failure index, \(k_\text{r}\) is a reference permeability, \(b\) is a fitting parameter. \(k_\text{r}\) and \(b\) can be calibrated by experimental data.
The failure index \(f\) is calculated from the Mohr Coulomb failure criterion comparing an acting shear stress for the shear dominated failure. The tensile failure is governed by an input parameter of tensile_strength_parameter .
The Mohr Coulomb failure criterion [25] takes the form
\[\tau(\sigma)=c-\sigma \mathrm{tan} \phi\]
with \(\tau\) the shear stress, \(c\) the cohesion, \(\sigma\) the normal stress, and \(\phi\) the internal friction angle.
The failure index of the Mohr Coulomb model is calculated by
\[ f_{MC}=\frac{|\tau_m| }{\cos(\phi)\tau(\sigma_m)} \]
with \(\tau_m=(\sigma_3-\sigma_1)/2\) and \(\sigma_m=(\sigma_1+\sigma_3)/2\), where \(\sigma_1\) and \(\sigma_3\) are the minimum and maximum shear stress, respectively.
The tensile failure index is calculated by
\[ f_{t} = \sigma_m / \sigma^t_{max} \]
with, \(0 < \sigma^t_{max} < c \tan(\phi) \), a parameter of tensile strength for the cutting of the apex of the Mohr Coulomb model.
The tensile stress status is determined by a condition of \(\sigma_m> \sigma^t_{max}\). The failure index is then calculated by
\[ f = \begin{cases} f_{MC}, & \sigma_{m} \leq \sigma^t_{max}\\ max(f_{MC}, f_t), & \sigma_{m} > \sigma^t_{max}\\ \end{cases} \]
The computed permeability components are restricted with an upper bound, i.e. \(\mathbf{k}:=k_{ij} < k_{max}\).
If \(\mathbf{k}_0\) is orthogonal, i.e input two or three numbers for its diagonal entries, a coordinate system rotation of \(\mathbf{k}\) is possible if it is needed.
Note: the conventional mechanics notations are used, which mean that tensile stress is positive.
Definition at line 91 of file PermeabilityMohrCoulombFailureIndexModel.h.
#include <PermeabilityMohrCoulombFailureIndexModel.h>
Public Member Functions | |
PermeabilityMohrCoulombFailureIndexModel (std::string name, ParameterLib::Parameter< double > const &k0, double const kr, double const b, double const c, double const phi, double const k_max, double const t_sigma_max, ParameterLib::CoordinateSystem const *const local_coordinate_system) | |
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 | |
ParameterLib::Parameter< double > const & | k0_ |
double const | kr_ |
Reference permeability. | |
double const | b_ |
Fitting parameter. | |
double const | c_ |
Cohesion. | |
double const | phi_ |
Angle of internal friction. | |
double const | k_max_ |
Maximum permeability. | |
double const | t_sigma_max_ |
Tensile strength parameter. | |
ParameterLib::CoordinateSystem const *const | local_coordinate_system_ |
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_ |
MaterialPropertyLib::PermeabilityMohrCoulombFailureIndexModel< DisplacementDim >::PermeabilityMohrCoulombFailureIndexModel | ( | std::string | name, |
ParameterLib::Parameter< double > const & | k0, | ||
double const | kr, | ||
double const | b, | ||
double const | c, | ||
double const | phi, | ||
double const | k_max, | ||
double const | t_sigma_max, | ||
ParameterLib::CoordinateSystem const *const | local_coordinate_system ) |
Definition at line 31 of file PermeabilityMohrCoulombFailureIndexModel.cpp.
References MaterialPropertyLib::PermeabilityMohrCoulombFailureIndexModel< DisplacementDim >::c_, MaterialPropertyLib::name, MaterialPropertyLib::Property::name_, OGS_FATAL, MaterialPropertyLib::PermeabilityMohrCoulombFailureIndexModel< DisplacementDim >::phi_, and MaterialPropertyLib::PermeabilityMohrCoulombFailureIndexModel< DisplacementDim >::t_sigma_max_.
|
overridevirtual |
Reimplemented from MaterialPropertyLib::Property.
Definition at line 61 of file PermeabilityMohrCoulombFailureIndexModel.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 152 of file PermeabilityMohrCoulombFailureIndexModel.cpp.
References MaterialPropertyLib::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 74 of file PermeabilityMohrCoulombFailureIndexModel.cpp.
References MaterialPropertyLib::formEigenTensor< 3 >(), MaterialPropertyLib::fromVector(), and MaterialPropertyLib::VariableArray::total_stress.
|
private |
Fitting parameter.
Definition at line 117 of file PermeabilityMohrCoulombFailureIndexModel.h.
|
private |
Cohesion.
Definition at line 119 of file PermeabilityMohrCoulombFailureIndexModel.h.
Referenced by MaterialPropertyLib::PermeabilityMohrCoulombFailureIndexModel< DisplacementDim >::PermeabilityMohrCoulombFailureIndexModel().
|
private |
Intrinsic permeability for undamaged material. It can be a scalar or tensor for anisotropic material.
Definition at line 113 of file PermeabilityMohrCoulombFailureIndexModel.h.
|
private |
Maximum permeability.
Definition at line 124 of file PermeabilityMohrCoulombFailureIndexModel.h.
|
private |
Reference permeability.
Definition at line 115 of file PermeabilityMohrCoulombFailureIndexModel.h.
|
private |
Definition at line 128 of file PermeabilityMohrCoulombFailureIndexModel.h.
|
private |
Angle of internal friction.
Definition at line 121 of file PermeabilityMohrCoulombFailureIndexModel.h.
Referenced by MaterialPropertyLib::PermeabilityMohrCoulombFailureIndexModel< DisplacementDim >::PermeabilityMohrCoulombFailureIndexModel().
|
private |
Tensile strength parameter.
Definition at line 127 of file PermeabilityMohrCoulombFailureIndexModel.h.
Referenced by MaterialPropertyLib::PermeabilityMohrCoulombFailureIndexModel< DisplacementDim >::PermeabilityMohrCoulombFailureIndexModel().