KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
linear_elastic_plane_strain_2D_law.h
Go to the documentation of this file.
1 // KRATOS___
2 // // ) )
3 // // ___ ___
4 // // ____ //___) ) // ) )
5 // // / / // // / /
6 // ((____/ / ((____ ((___/ / MECHANICS
7 //
8 // License: geo_mechanics_application/license.txt
9 //
10 // Main authors: Vahid Galavi
11 //
12 
13 #pragma once
14 
15 // Project includes
17 
18 namespace Kratos
19 {
20 
28 class KRATOS_API(GEO_MECHANICS_APPLICATION) GeoLinearElasticPlaneStrain2DLaw
29  : public LinearPlaneStrainK0Law
30 {
31 public:
34 
36  using SizeType = std::size_t;
37 
39  static constexpr SizeType Dimension = N_DIM_2D;
40 
42  // for the time being
43  static constexpr SizeType VoigtSize = VOIGT_SIZE_2D_PLANE_STRAIN;
44 
47 
51  ConstitutiveLaw::Pointer Clone() const override;
52 
53  bool RequiresInitializeMaterialResponse() override;
54  void InitializeMaterialResponseCauchy(ConstitutiveLaw::Parameters& rValues) override;
55 
56  bool RequiresFinalizeMaterialResponse() override;
57  void FinalizeMaterialResponseCauchy(ConstitutiveLaw::Parameters & rValues) override;
58  void FinalizeMaterialResponsePK2(ConstitutiveLaw::Parameters& rValues) override;
59 
64  void GetLawFeatures(Features& rFeatures) override;
65 
70  SizeType WorkingSpaceDimension() override;
71 
76  SizeType GetStrainSize() const override;
77 
78  bool IsIncremental() override;
79 
86  bool& GetValue(const Variable<bool>& rThisVariable, bool& rValue) override;
88 
89 protected:
95  void CalculateElasticMatrix(Matrix& C, ConstitutiveLaw::Parameters& rValues) override;
96 
103  void CalculatePK2Stress(const Vector& rStrainVector,
104  Vector& rStressVector,
105  ConstitutiveLaw::Parameters& rValues) override;
106 
107 
108 
110 
111 private:
112  Vector mStressVector = ZeroVector(VoigtSize);
113  Vector mStressVectorFinalized = ZeroVector(VoigtSize);
114  Vector mDeltaStrainVector = ZeroVector(VoigtSize);
115  Vector mStrainVectorFinalized = ZeroVector(VoigtSize);
116  bool mIsModelInitialized = false;
117 
118  friend class Serializer;
119  void save(Serializer& rSerializer) const override;
120  void load(Serializer& rSerializer) override;
121 }; // Class GeoLinearElasticPlaneStrain2DLaw
122 
123 }
virtual bool & GetValue(const Variable< bool > &rThisVariable, bool &rValue)
Returns the value of a specified variable (boolean)
Definition: constitutive_law.cpp:201
std::size_t SizeType
Definition: constitutive_law.h:82
This class defines a small deformation linear elastic constitutive model for plane strain cases.
Definition: linear_elastic_plane_strain_2D_law.h:30
KRATOS_CLASS_POINTER_DEFINITION(GeoLinearElasticPlaneStrain2DLaw)
Counted pointer of LinearPlaneStrainK0Law.
This class defines a small deformation linear elastic constitutive model for plane strain cases.
Definition: linear_elastic_plane_strain_K0_law.h:49
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
Parameters GetValue(Parameters &rParameters, const std::string &rEntry)
Definition: add_kratos_parameters_to_python.cpp:53
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
KratosZeroVector< double > ZeroVector
Definition: amatrix_interface.h:561
constexpr SizeType VOIGT_SIZE_2D_PLANE_STRAIN
Definition: geo_mechanics_application_constants.h:43
constexpr SizeType N_DIM_2D
Definition: geo_mechanics_application_constants.h:26
int C
Definition: generate_hyper_elastic_simo_taylor_neo_hookean.py:27
def load(f)
Definition: ode_solve.py:307
Definition: constitutive_law.h:137
Definition: constitutive_law.h:189