OGS
MaterialPropertyLib::PermeabilityMohrCoulombFailureIndexModel< DisplacementDim > Class Template Referencefinal

Detailed Description

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

A failure index dependent permeability model [44].

\[ \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 [24] 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 83 of file PermeabilityMohrCoulombFailureIndexModel.h.

#include <PermeabilityMohrCoulombFailureIndexModel.h>

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

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>
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

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_

Constructor & Destructor Documentation

◆ PermeabilityMohrCoulombFailureIndexModel()

template<int DisplacementDim>
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 23 of file PermeabilityMohrCoulombFailureIndexModel.cpp.

29 : k0_(k0),
30 kr_(kr),
31 b_(b),
32 c_(c),
37{
38 const double t_sigma_upper = c_ / std::tan(phi_);
42 {
44 "Tensile strength parameter of {:e} is out of the range (0, "
45 "c/tan(phi)) = (0, {:e})",
47 }
48
50}
#define OGS_FATAL(...)
Definition Error.h:19
A failure index dependent permeability model wangetalperm2020.

References b_, MaterialPropertyLib::c, c_, k0_, k_max_, kr_, local_coordinate_system_, MaterialPropertyLib::name, MaterialPropertyLib::Property::name_, OGS_FATAL, phi_, and t_sigma_max_.

Member Function Documentation

◆ checkScale()

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

Reimplemented from MaterialPropertyLib::Property.

Definition at line 53 of file PermeabilityMohrCoulombFailureIndexModel.cpp.

55{
57 {
59 "The property 'PermeabilityMohrCoulombFailureIndexModel' is "
60 "implemented on the 'medium' scale only.");
61 }
62}
std::variant< Medium *, Phase *, Component * > scale_

References OGS_FATAL, and MaterialPropertyLib::Property::scale_.

◆ dValue()

template<int DisplacementDim>
PropertyDataType MaterialPropertyLib::PermeabilityMohrCoulombFailureIndexModel< 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 144 of file PermeabilityMohrCoulombFailureIndexModel.cpp.

148{
150 {
151 return 0.;
152 }
153
154 OGS_FATAL(
155 "The derivative of the intrinsic permeability k(sigma, ...) with "
156 "respect to stress tensor (sigma) is not implemented because that "
157 "dk/du is normally omitted.");
158}

References MaterialPropertyLib::mechanical_strain, and OGS_FATAL.

◆ value()

template<int DisplacementDim>
PropertyDataType MaterialPropertyLib::PermeabilityMohrCoulombFailureIndexModel< 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 66 of file PermeabilityMohrCoulombFailureIndexModel.cpp.

70{
71 auto const& stress_vector =
73
74 auto const& stress_tensor =
76
79
80 // Principle stress
81 auto const sigma = eigenvalue_solver.eigenvalues();
82
83 auto k_data = k0_(t, pos);
84
85 double const max_sigma = std::max(std::fabs(sigma[0]), std::fabs(sigma[2]));
86
88 {
89 return fromVector(k_data);
90 }
91
92 double const sigma_m = 0.5 * (sigma[2] + sigma[0]);
93
94 double const tau_m = 0.5 * std::fabs(sigma[2] - sigma[0]);
95 double f = 0.0;
97 {
98 // tensile failure criterion
100
101 double const tau_tt =
103
104 f = std::max(f, tau_m / tau_tt);
105 }
106 else
107 {
108 // Mohr Coulomb failure criterion
109 f = tau_m / (c_ * std::cos(phi_) - sigma_m * std::sin(phi_));
110 }
111
112 if (f >= 1.0)
113 {
114 const double exp_value = std::exp(b_ * f);
116 [&](double const k_i)
117 { return std::min(k_i + kr_ * exp_value, k_max_); });
118 }
119
120 // Local coordinate transformation is only applied for the case that the
121 // initial intrinsic permeability is given with orthotropic assumption.
123 {
128
129 for (int i = 0; i < DisplacementDim; ++i)
130 {
132 ei_otimes_ei = e.col(i) * e.col(i).transpose();
133
134 k += k_data[i] * ei_otimes_ei;
135 }
136 return k;
137 }
138
139 return fromVector(k_data);
140}
template Eigen::Matrix< double, 3, 3 > formEigenTensor< 3 >(MaterialPropertyLib::PropertyDataType const &values)
PropertyDataType fromVector(std::vector< double > const &values)

References b_, c_, MaterialPropertyLib::formEigenTensor< 3 >(), MaterialPropertyLib::fromVector(), k0_, k_max_, kr_, local_coordinate_system_, phi_, t_sigma_max_, and MaterialPropertyLib::VariableArray::total_stress.

Member Data Documentation

◆ b_

template<int DisplacementDim>
double const MaterialPropertyLib::PermeabilityMohrCoulombFailureIndexModel< DisplacementDim >::b_
private

Fitting parameter.

Definition at line 109 of file PermeabilityMohrCoulombFailureIndexModel.h.

Referenced by PermeabilityMohrCoulombFailureIndexModel(), and value().

◆ c_

template<int DisplacementDim>
double const MaterialPropertyLib::PermeabilityMohrCoulombFailureIndexModel< DisplacementDim >::c_
private

Cohesion.

Definition at line 111 of file PermeabilityMohrCoulombFailureIndexModel.h.

Referenced by PermeabilityMohrCoulombFailureIndexModel(), and value().

◆ k0_

template<int DisplacementDim>
ParameterLib::Parameter<double> const& MaterialPropertyLib::PermeabilityMohrCoulombFailureIndexModel< DisplacementDim >::k0_
private

Intrinsic permeability for undamaged material. It can be a scalar or tensor for anisotropic material.

Definition at line 105 of file PermeabilityMohrCoulombFailureIndexModel.h.

Referenced by PermeabilityMohrCoulombFailureIndexModel(), and value().

◆ k_max_

template<int DisplacementDim>
double const MaterialPropertyLib::PermeabilityMohrCoulombFailureIndexModel< DisplacementDim >::k_max_
private

Maximum permeability.

Definition at line 116 of file PermeabilityMohrCoulombFailureIndexModel.h.

Referenced by PermeabilityMohrCoulombFailureIndexModel(), and value().

◆ kr_

template<int DisplacementDim>
double const MaterialPropertyLib::PermeabilityMohrCoulombFailureIndexModel< DisplacementDim >::kr_
private

Reference permeability.

Definition at line 107 of file PermeabilityMohrCoulombFailureIndexModel.h.

Referenced by PermeabilityMohrCoulombFailureIndexModel(), and value().

◆ local_coordinate_system_

template<int DisplacementDim>
ParameterLib::CoordinateSystem const* const MaterialPropertyLib::PermeabilityMohrCoulombFailureIndexModel< DisplacementDim >::local_coordinate_system_
private

◆ phi_

template<int DisplacementDim>
double const MaterialPropertyLib::PermeabilityMohrCoulombFailureIndexModel< DisplacementDim >::phi_
private

Angle of internal friction.

Definition at line 113 of file PermeabilityMohrCoulombFailureIndexModel.h.

Referenced by PermeabilityMohrCoulombFailureIndexModel(), and value().

◆ t_sigma_max_

template<int DisplacementDim>
double const MaterialPropertyLib::PermeabilityMohrCoulombFailureIndexModel< DisplacementDim >::t_sigma_max_
private

Tensile strength parameter.

Definition at line 119 of file PermeabilityMohrCoulombFailureIndexModel.h.

Referenced by PermeabilityMohrCoulombFailureIndexModel(), and value().


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