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_plane_stress.h
Go to the documentation of this file.
1 // KRATOS ___| | | |
2 // \___ \ __| __| | | __| __| | | __| _` | |
3 // | | | | | ( | | | | ( | |
4 // _____/ \__|_| \__,_|\___|\__|\__,_|_| \__,_|_| MECHANICS
5 //
6 // License: BSD License
7 // license: StructuralMechanicsApplication/license.txt
8 //
9 // Main authors: Riccardo Rossi
10 //
11 
12 #pragma once
13 
14 // System includes
15 
16 // External includes
17 
18 // Project includes
19 #include "custom_constitutive/elastic_isotropic_3d.h"
20 
21 namespace Kratos
22 {
25 
29 
33 
37 
41 
48 class KRATOS_API(STRUCTURAL_MECHANICS_APPLICATION) LinearPlaneStress
49  : public ElasticIsotropic3D
50 {
51 public:
54 
57 
60 
63 
64  // Adding the respective using to avoid overload conflicts
65  using BaseType::Has;
66  using BaseType::GetValue;
67 
69  typedef std::size_t SizeType;
70 
72  static constexpr SizeType Dimension = 2;
73 
75  static constexpr SizeType VoigtSize = 3;
76 
79 
82 
87 
88  ConstitutiveLaw::Pointer Clone() const override;
89 
93  LinearPlaneStress (const LinearPlaneStress& rOther);
94 
95 
99  ~LinearPlaneStress() override;
100 
104 
108 
113  void GetLawFeatures(Features& rFeatures) override;
114 
119  {
120  return Dimension;
121  };
122 
126  SizeType GetStrainSize() const override
127  {
128  return VoigtSize;
129  }
130 
134 
138 
142 
146 
153  bool& GetValue(const Variable<bool>& rThisVariable, bool& rValue) override;
154 
156 
157 protected:
158 
161 
165 
169 
173 
179  void CalculateElasticMatrix(VoigtSizeMatrixType& C, ConstitutiveLaw::Parameters& rValues) override;
180 
187  void CalculatePK2Stress(
188  const ConstitutiveLaw::StrainVectorType& rStrainVector,
189  ConstitutiveLaw::StressVectorType& rStressVector,
191  ) override;
192 
198  void CalculateCauchyGreenStrain(
200  ConstitutiveLaw::StrainVectorType& rStrainVector
201  ) override;
202 
204 
205 private:
206 
209 
213 
217 
222 
227 
231  friend class Serializer;
232 
233  void save(Serializer& rSerializer) const override
234  {
236  }
237 
238  void load(Serializer& rSerializer) override
239  {
241  }
242 
243 
244 }; // Class LinearPlaneStress
245 } // namespace Kratos.
Definition: constitutive_law.h:47
Definition: elastic_isotropic_3d.h:53
This class defines a small deformation linear elastic constitutive model for plane stress cases.
Definition: linear_plane_stress.h:50
std::size_t SizeType
The size type definition.
Definition: linear_plane_stress.h:69
KRATOS_CLASS_POINTER_DEFINITION(LinearPlaneStress)
Counted pointer of LinearPlaneStress.
SizeType WorkingSpaceDimension() override
Definition: linear_plane_stress.h:118
SizeType GetStrainSize() const override
Definition: linear_plane_stress.h:126
ProcessInfo ProcessInfoType
The process info definition.
Definition: linear_plane_stress.h:56
ElasticIsotropic3D BaseType
The base class ElasticIsotropic3D type definition.
Definition: linear_plane_stress.h:62
ConstitutiveLaw CLBaseType
The base class ConstitutiveLaw type definition.
Definition: linear_plane_stress.h:59
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
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
#define KRATOS_SERIALIZE_SAVE_BASE_CLASS(Serializer, BaseType)
Definition: define.h:812
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
bool Has(const std::string &ModelerName)
Checks if the modeler is registered.
Definition: modeler_factory.cpp:24
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
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