OGS
MaterialPropertyLib::EmbeddedFracturePermeability< DisplacementDim > Class Template Referencefinal

Detailed Description

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

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 [23] .

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 49 of file EmbeddedFracturePermeability.h.

#include <EmbeddedFracturePermeability.h>

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

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)
 
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 primary_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. More...
 
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
 

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
 

Additional Inherited Members

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

Member Typedef Documentation

◆ SymmetricTensor

template<int DisplacementDim>
using MaterialPropertyLib::EmbeddedFracturePermeability< DisplacementDim >::SymmetricTensor = Eigen::Matrix<double, KelvinVectorSize, 1>

Definition at line 72 of file EmbeddedFracturePermeability.h.

Constructor & Destructor Documentation

◆ EmbeddedFracturePermeability()

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

Definition at line 20 of file EmbeddedFracturePermeability.cpp.

28  : _n(fracture_normal),
29  _n_const(fracture_normal_is_constant),
30  _k(intrinsic_permeability),
31  _b0(initial_aperture),
32  _a(mean_fracture_distance),
33  _e0(threshold_strain)
34 {
35  name_ = std::move(name);
36 }

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

Member Function Documentation

◆ checkScale()

template<int DisplacementDim>
void MaterialPropertyLib::EmbeddedFracturePermeability< DisplacementDim >::checkScale
overridevirtual

Reimplemented from MaterialPropertyLib::Property.

Definition at line 39 of file EmbeddedFracturePermeability.cpp.

40 {
41  if (!std::holds_alternative<Medium*>(scale_))
42  {
43  OGS_FATAL(
44  "The property 'EmbeddedFracturePermeability' is "
45  "implemented on the 'medium' scale only.");
46  }
47 }
#define OGS_FATAL(...)
Definition: Error.h:26
std::variant< Medium *, Phase *, Component * > scale_
Definition: Property.h:287

References OGS_FATAL.

◆ dValue()

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

86 {
87  if (primary_variable != Variable::mechanical_strain)
88  {
89  OGS_FATAL(
90  "EmbeddedFracturePermeability::dValue is implemented for "
91  "derivatives with respect to strain only.");
92  }
93 
94  Eigen::Matrix<double, 3, 1> const n = [&]
95  {
96  if (_n_const)
97  {
98  return _n;
99  }
100  auto const sigma = formEigenTensor<3>(std::get<SymmetricTensor>(
101  variable_array[static_cast<int>(Variable::total_stress)]));
102  Eigen::SelfAdjointEigenSolver<Eigen::Matrix<double, 3, 3>> e_s(sigma);
103  return (Eigen::Matrix<double, 3, 1>)e_s.eigenvectors().col(2);
104  }();
105 
108  variable_array[static_cast<int>(Variable::mechanical_strain)]));
109  double const e_n = (eps * n).dot(n.transpose());
110  double const H_de = (e_n > _e0) ? 1.0 : 0.0;
111  double const b_f = _b0 + H_de * _a * (e_n - _e0);
112 
113  Eigen::Matrix3d const M = n * n.transpose();
114  return (H_de * (b_f * b_f / 4 - _k) * (Eigen::Matrix3d::Identity() - M) * M)
115  .template topLeftCorner<DisplacementDim, DisplacementDim>()
116  .eval();
117 }
template Eigen::Matrix< double, 3, 3 > formEigenTensor< 3 >(MaterialPropertyLib::PropertyDataType const &values)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Definition: KelvinVector.h:48
Eigen::Matrix< double, 3, 3 > kelvinVectorToTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)

References MaterialPropertyLib::formEigenTensor< 3 >(), MathLib::KelvinVector::kelvinVectorToTensor(), MaterialPropertyLib::mechanical_strain, OGS_FATAL, and MaterialPropertyLib::total_stress.

◆ value()

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

54 {
55  Eigen::Matrix<double, 3, 1> const n = [&]
56  {
57  if (_n_const)
58  {
59  return _n;
60  }
61  auto const sigma = formEigenTensor<3>(std::get<SymmetricTensor>(
62  variable_array[static_cast<int>(Variable::total_stress)]));
63  Eigen::SelfAdjointEigenSolver<Eigen::Matrix<double, 3, 3>> e_s(sigma);
64  return (Eigen::Matrix<double, 3, 1>)e_s.eigenvectors().col(2);
65  }();
66 
69  variable_array[static_cast<int>(Variable::mechanical_strain)]));
70  double const e_n = (eps * n).dot(n.transpose());
71  double const H_de = (e_n > _e0) ? 1.0 : 0.0;
72  double const b_f = _b0 + H_de * _a * (e_n - _e0);
73  double const coeff = H_de * (b_f / _a) * ((b_f * b_f / 12.0) - _k);
74 
75  Eigen::Matrix3d I = Eigen::Matrix3d::Identity();
76  return (coeff * (I - n * n.transpose()) + _k * I)
77  .template topLeftCorner<DisplacementDim, DisplacementDim>()
78  .eval();
79 }

References MaterialPropertyLib::formEigenTensor< 3 >(), MathLib::KelvinVector::kelvinVectorToTensor(), MaterialPropertyLib::mechanical_strain, and MaterialPropertyLib::total_stress.

Member Data Documentation

◆ _a

template<int DisplacementDim>
double const MaterialPropertyLib::EmbeddedFracturePermeability< DisplacementDim >::_a
private

Definition at line 57 of file EmbeddedFracturePermeability.h.

◆ _b0

template<int DisplacementDim>
double const MaterialPropertyLib::EmbeddedFracturePermeability< DisplacementDim >::_b0
private

Definition at line 56 of file EmbeddedFracturePermeability.h.

◆ _e0

template<int DisplacementDim>
double const MaterialPropertyLib::EmbeddedFracturePermeability< DisplacementDim >::_e0
private

Definition at line 58 of file EmbeddedFracturePermeability.h.

◆ _k

template<int DisplacementDim>
double const MaterialPropertyLib::EmbeddedFracturePermeability< DisplacementDim >::_k
private

Definition at line 55 of file EmbeddedFracturePermeability.h.

◆ _medium

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

Definition at line 52 of file EmbeddedFracturePermeability.h.

◆ _n

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

Definition at line 53 of file EmbeddedFracturePermeability.h.

◆ _n_const

template<int DisplacementDim>
bool const MaterialPropertyLib::EmbeddedFracturePermeability< DisplacementDim >::_n_const
private

Definition at line 54 of file EmbeddedFracturePermeability.h.

◆ KelvinVectorSize

template<int DisplacementDim>
int const MaterialPropertyLib::EmbeddedFracturePermeability< DisplacementDim >::KelvinVectorSize
static
Initial value:
=
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
Definition: KelvinVector.h:23

Definition at line 70 of file EmbeddedFracturePermeability.h.


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