OGS
MaterialLib::Fracture::LinearElasticIsotropic< DisplacementDim > Class Template Referencefinal

Detailed Description

template<int DisplacementDim>
class MaterialLib::Fracture::LinearElasticIsotropic< DisplacementDim >

Definition at line 24 of file LinearElasticIsotropic.h.

#include <LinearElasticIsotropic.h>

Inheritance diagram for MaterialLib::Fracture::LinearElasticIsotropic< DisplacementDim >:
[legend]
Collaboration diagram for MaterialLib::Fracture::LinearElasticIsotropic< DisplacementDim >:
[legend]

Classes

struct  MaterialProperties
 Variables specific to the material model. More...
struct  MaterialStateVariables

Public Member Functions

std::unique_ptr< typename FractureModelBase< DisplacementDim >::MaterialStateVariables > createMaterialStateVariables () override
 LinearElasticIsotropic (double const penalty_aperture_cutoff, bool const tension_cutoff, MaterialProperties material_properties)
void computeConstitutiveRelation (double const t, ParameterLib::SpatialPosition const &x, double const aperture0, Eigen::Ref< Eigen::VectorXd const > sigma0, Eigen::Ref< Eigen::VectorXd const > w_prev, Eigen::Ref< Eigen::VectorXd const > w, Eigen::Ref< Eigen::VectorXd const > sigma_prev, Eigen::Ref< Eigen::VectorXd > sigma, Eigen::Ref< Eigen::MatrixXd > C, typename FractureModelBase< DisplacementDim >::MaterialStateVariables &material_state_variables) override
Public Member Functions inherited from MaterialLib::Fracture::FractureModelBase< DisplacementDim >
virtual ~FractureModelBase ()=default
virtual void computeConstitutiveRelation (double const t, ParameterLib::SpatialPosition const &x, double const aperture0, Eigen::Ref< Eigen::VectorXd const > sigma0, Eigen::Ref< Eigen::VectorXd const > w_prev, Eigen::Ref< Eigen::VectorXd const > w, Eigen::Ref< Eigen::VectorXd const > sigma_prev, Eigen::Ref< Eigen::VectorXd > sigma, Eigen::Ref< Eigen::MatrixXd > C, MaterialStateVariables &material_state_variables)=0

Private Attributes

double const _penalty_aperture_cutoff
bool const _tension_cutoff
MaterialProperties _mp

Constructor & Destructor Documentation

◆ LinearElasticIsotropic()

template<int DisplacementDim>
MaterialLib::Fracture::LinearElasticIsotropic< DisplacementDim >::LinearElasticIsotropic ( double const penalty_aperture_cutoff,
bool const tension_cutoff,
MaterialProperties material_properties )
inlineexplicit

Member Function Documentation

◆ computeConstitutiveRelation()

template<int DisplacementDim>
void MaterialLib::Fracture::LinearElasticIsotropic< DisplacementDim >::computeConstitutiveRelation ( double const t,
ParameterLib::SpatialPosition const & x,
double const aperture0,
Eigen::Ref< Eigen::VectorXd const > sigma0,
Eigen::Ref< Eigen::VectorXd const > w_prev,
Eigen::Ref< Eigen::VectorXd const > w,
Eigen::Ref< Eigen::VectorXd const > sigma_prev,
Eigen::Ref< Eigen::VectorXd > sigma,
Eigen::Ref< Eigen::MatrixXd > C,
typename FractureModelBase< DisplacementDim >::MaterialStateVariables & material_state_variables )
override

Computation of the constitutive relation for the linear elastic model.

Parameters
tcurrent time
xcurrent position in space
aperture0initial fracture's aperture
sigma0initial stress
w_prevfracture displacement at previous time step
wfracture displacement at current time step
sigma_prevstress at previous time step
sigmastress at current time step
Ctangent matrix for stress and fracture displacements
material_state_variablesmaterial state variables

Definition at line 19 of file LinearElasticIsotropic.cpp.

37{
39
40 const int index_ns = DisplacementDim - 1;
41 C.setZero();
42 for (int i = 0; i < index_ns; i++)
43 {
44 C(i, i) = _mp.shear_stiffness(t, x)[0];
45 }
46
47 sigma.noalias() = C * w;
48
49 double const aperture = w[index_ns] + aperture0;
50
51 sigma.coeffRef(index_ns) =
52 _mp.normal_stiffness(t, x)[0] * w[index_ns] *
54
56 _mp.normal_stiffness(t, x)[0] *
58
59 sigma.noalias() += sigma0;
60
61 // correction for an opening fracture
62 if (_tension_cutoff && sigma[index_ns] > 0)
63 {
64 C.setZero();
65 sigma.setZero();
66 material_state_variables.setTensileStress(true);
67 }
68}
double logPenalty(double const aperture0, double const aperture, double const aperture_cutoff)
Definition LogPenalty.h:49
double logPenaltyDerivative(double const aperture0, double const aperture, double const aperture_cutoff)
Definition LogPenalty.h:23

References _mp, _penalty_aperture_cutoff, _tension_cutoff, MaterialLib::Fracture::logPenalty(), MaterialLib::Fracture::logPenaltyDerivative(), MaterialLib::Fracture::FractureModelBase< DisplacementDim >::MaterialStateVariables::reset(), and MaterialLib::Fracture::FractureModelBase< DisplacementDim >::MaterialStateVariables::setTensileStress().

◆ createMaterialStateVariables()

template<int DisplacementDim>
std::unique_ptr< typename FractureModelBase< DisplacementDim >::MaterialStateVariables > MaterialLib::Fracture::LinearElasticIsotropic< DisplacementDim >::createMaterialStateVariables ( )
inlineoverridevirtual

Member Data Documentation

◆ _mp

template<int DisplacementDim>
MaterialProperties MaterialLib::Fracture::LinearElasticIsotropic< DisplacementDim >::_mp
private

◆ _penalty_aperture_cutoff

template<int DisplacementDim>
double const MaterialLib::Fracture::LinearElasticIsotropic< DisplacementDim >::_penalty_aperture_cutoff
private

Compressive normal displacements above this value will not enter the computation of the normal stiffness modulus of the fracture.

Note
Setting this to the initial aperture value allows negative apertures.

Definition at line 107 of file LinearElasticIsotropic.h.

Referenced by LinearElasticIsotropic(), and computeConstitutiveRelation().

◆ _tension_cutoff

template<int DisplacementDim>
bool const MaterialLib::Fracture::LinearElasticIsotropic< DisplacementDim >::_tension_cutoff
private

If set no resistance to open the fracture over the initial aperture is opposed.

Definition at line 111 of file LinearElasticIsotropic.h.

Referenced by LinearElasticIsotropic(), and computeConstitutiveRelation().


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