OGS
CohesiveZoneModeI.h
Go to the documentation of this file.
1
10#pragma once
11
12#include <Eigen/Core>
13#include <utility>
14
15#include "FractureModelBase.h"
17
18namespace MaterialLib
19{
20namespace Fracture
21{
22namespace CohesiveZoneModeI
23{
26{
29
30 MaterialPropertiesParameters(P const& normal_stiffness_,
31 P const& shear_stiffness_,
32 P const& fracture_toughness_,
33 P const& peak_normal_traction_)
34 : normal_stiffness(normal_stiffness_),
35 shear_stiffness(shear_stiffness_),
36 fracture_toughness(fracture_toughness_),
37 peak_normal_traction(peak_normal_traction_)
38 {
39 }
40
42 double fracture_opening_at_peak_traction(double const t, X const& x) const
43 {
44 return peak_normal_traction(t, x)[0] / normal_stiffness(t, x)[0];
45 }
46
49 X const& x) const
50 {
51 if (peak_normal_traction(t, x)[0] == 0.0)
52 {
53 return 0.0;
54 }
55
56 return 2 * fracture_toughness(t, x)[0] / peak_normal_traction(t, x)[0];
57 }
58
63
80};
81
85{
88 : Kn(mp.normal_stiffness(t, x)[0]),
89 Ks(mp.shear_stiffness(t, x)[0]),
90 Gc(mp.fracture_toughness(t, x)[0]),
91 t_np(mp.peak_normal_traction(t, x)[0]),
92 w_np(mp.fracture_opening_at_peak_traction(t, x)),
93 w_nf(mp.fracture_opening_at_residual_traction(t, x))
94 {
95 }
96
97 double const Kn;
98 double const Ks;
99 double const Gc;
100 double const t_np;
101 double const w_np;
102 double const w_nf;
103};
104
105template <int DisplacementDim>
107 : public FractureModelBase<DisplacementDim>::MaterialStateVariables
108{
110
111 void pushBackState() override { damage_prev = damage; }
112
113 double damage;
114
115 // Initial values from previous timestep
116 double damage_prev;
117};
118
124template <int DisplacementDim>
125class CohesiveZoneModeI final : public FractureModelBase<DisplacementDim>
126{
127public:
128 std::unique_ptr<
131 {
132 return std::make_unique<StateVariables<DisplacementDim>>();
133 }
134
135public:
136public:
137 explicit CohesiveZoneModeI(double const penalty_aperture_cutoff,
138 bool const tension_cutoff,
139 MaterialPropertiesParameters material_properties)
140 : _penalty_aperture_cutoff(penalty_aperture_cutoff),
141 _tension_cutoff(tension_cutoff),
142 _mp(std::move(material_properties))
143 {
144 }
145
162 double const t,
164 double const aperture0,
165 Eigen::Ref<Eigen::VectorXd const>
166 sigma0,
167 Eigen::Ref<Eigen::VectorXd const>
168 w_prev,
169 Eigen::Ref<Eigen::VectorXd const>
170 w,
171 Eigen::Ref<Eigen::VectorXd const>
172 sigma_prev,
173 Eigen::Ref<Eigen::VectorXd>
174 sigma,
175 Eigen::Ref<Eigen::MatrixXd>
176 C,
178 material_state_variables) override;
179
181 double const t, ParameterLib::SpatialPosition const& x) const
182 {
183 return MaterialProperties(t, x, _mp);
184 }
185
186private:
192
195 bool const _tension_cutoff;
196
198};
199
200} // namespace CohesiveZoneModeI
201} // namespace Fracture
202} // namespace MaterialLib
203
204namespace MaterialLib
205{
206namespace Fracture
207{
208namespace CohesiveZoneModeI
209{
210extern template class CohesiveZoneModeI<2>;
211extern template class CohesiveZoneModeI<3>;
212} // namespace CohesiveZoneModeI
213} // namespace Fracture
214} // namespace MaterialLib
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
CohesiveZoneModeI(double const penalty_aperture_cutoff, bool const tension_cutoff, MaterialPropertiesParameters material_properties)
MaterialProperties evaluatedMaterialProperties(double const t, ParameterLib::SpatialPosition const &x) const
std::unique_ptr< typename FractureModelBase< DisplacementDim >::MaterialStateVariables > createMaterialStateVariables() override
P const & shear_stiffness
Shear stiffness given in units of stress per length.
double fracture_opening_at_peak_traction(double const t, X const &x) const
Assuming initially stress-free state.
double fracture_opening_at_residual_traction(double const t, X const &x) const
Assuming initially stress-free state.
P const & normal_stiffness
Normal stiffness given in units of stress per length.
MaterialPropertiesParameters(P const &normal_stiffness_, P const &shear_stiffness_, P const &fracture_toughness_, P const &peak_normal_traction_)
MaterialProperties(double const t, ParameterLib::SpatialPosition const &x, MaterialPropertiesParameters const &mp)