OGS
FractureModels/LinearElasticIsotropic.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
6#include <Eigen/Core>
7#include <utility>
8
9#include "FractureModelBase.h"
11
12namespace MaterialLib
13{
14namespace Fracture
15{
16
17template <int DisplacementDim>
18class LinearElasticIsotropic final : public FractureModelBase<DisplacementDim>
19{
20public:
23 {
26
27 MaterialProperties(P const& normal_stiffness_, P const& shear_stiffness_)
28 : normal_stiffness(normal_stiffness_), shear_stiffness(shear_stiffness_)
29 {
30 }
31
36 };
37
39 : public FractureModelBase<DisplacementDim>::MaterialStateVariables
40 {
41 void pushBackState() override {}
42 };
43
44 std::unique_ptr<
47 {
48 return std::unique_ptr<typename FractureModelBase<
49 DisplacementDim>::MaterialStateVariables>{
51 }
52
53public:
54 explicit LinearElasticIsotropic(double const penalty_aperture_cutoff,
55 bool const tension_cutoff,
56 MaterialProperties material_properties)
57 : _penalty_aperture_cutoff(penalty_aperture_cutoff),
58 _tension_cutoff(tension_cutoff),
59 _mp(std::move(material_properties))
60 {
61 }
62
78 double const t,
80 double const aperture0,
81 Eigen::Ref<Eigen::VectorXd const>
82 sigma0,
83 Eigen::Ref<Eigen::VectorXd const>
84 w_prev,
85 Eigen::Ref<Eigen::VectorXd const>
86 w,
87 Eigen::Ref<Eigen::VectorXd const>
88 sigma_prev,
89 Eigen::Ref<Eigen::VectorXd>
90 sigma,
91 Eigen::Ref<Eigen::MatrixXd>
92 C,
94 material_state_variables) override;
95
96private:
102
105 bool const _tension_cutoff;
106
108};
109
110} // namespace Fracture
111} // namespace MaterialLib
112
113namespace MaterialLib
114{
115namespace Fracture
116{
117extern template class LinearElasticIsotropic<2>;
118extern template class LinearElasticIsotropic<3>;
119} // namespace Fracture
120} // namespace MaterialLib
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
std::unique_ptr< typename FractureModelBase< DisplacementDim >::MaterialStateVariables > createMaterialStateVariables() override
P const & shear_stiffness
Shear stiffness given in units of stress per length.
P const & normal_stiffness
Normal stiffness given in units of stress per length.