25template <
int DisplacementDim>
35 Eigen::Matrix<double, DisplacementDim, 1>
w_p =
36 Eigen::Matrix<double, DisplacementDim, 1>::Zero();
39 Eigen::Matrix<double, DisplacementDim, 1>
w_p_prev;
44template <
int DisplacementDim>
52 return std::make_unique<StateVariables<DisplacementDim>>();
63 P const& normal_stiffness_,
P const& shear_stiffness_,
64 P const& friction_angle_,
P const& dilatancy_angle_,
91 double const penalty_aperture_cutoff,
92 bool const tension_cutoff,
97 _mp(std::move(material_properties))
118 double const aperture0,
119 Eigen::Ref<Eigen::VectorXd const>
121 Eigen::Ref<Eigen::VectorXd const>
123 Eigen::Ref<Eigen::VectorXd const>
125 Eigen::Ref<Eigen::VectorXd const>
127 Eigen::Ref<Eigen::VectorXd>
129 Eigen::Ref<Eigen::MatrixXd>
132 material_state_variables)
override;
160extern template class Coulomb<2>;
161extern template class Coulomb<3>;
bool const _tension_cutoff
NumLib::NewtonRaphsonSolverParameters const _nonlinear_solver_parameters
Coulomb(NumLib::NewtonRaphsonSolverParameters nonlinear_solver_parameters, double const penalty_aperture_cutoff, bool const tension_cutoff, MaterialProperties material_properties)
double const _penalty_aperture_cutoff
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 > Kep, typename FractureModelBase< DisplacementDim >::MaterialStateVariables &material_state_variables) override
std::unique_ptr< typename FractureModelBase< DisplacementDim >::MaterialStateVariables > createMaterialStateVariables() override
Variables specific to the material model.
P const & dilatancy_angle
P const & shear_stiffness
Shear stiffness given in units of stress per length.
P const & cohesion
Fracture cohesion in units of stress.
P const & normal_stiffness
Normal stiffness given in units of stress per length.
MaterialProperties(P const &normal_stiffness_, P const &shear_stiffness_, P const &friction_angle_, P const &dilatancy_angle_, P const &cohesion_)
void setInitialConditions()
Eigen::Matrix< double, DisplacementDim, 1 > w_p
void pushBackState() override
Eigen::Matrix< double, DisplacementDim, 1 > w_p_prev
EIGEN_MAKE_ALIGNED_OPERATOR_NEW