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.
axisym_elastic_isotropic.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: Vicente Mataix Ferrandiz
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 
42 class KRATOS_API(STRUCTURAL_MECHANICS_APPLICATION) AxisymElasticIsotropic : public ElasticIsotropic3D
43 {
44 public:
45 
48 
51  typedef std::size_t SizeType;
52 
57 
60 
65 
66  ConstitutiveLaw::Pointer Clone() const override;
67 
72 
76  ~AxisymElasticIsotropic() override;
77 
81 
85 
90  void GetLawFeatures(Features& rFeatures) override;
91 
95  SizeType GetStrainSize() const override
96  {
97  return 4;
98  }
99 
103 
107 
111 
115 
117 
118 protected:
119 
122 
126 
130 
134 
136 
137 private:
138 
144 
145 
149 
156  void CalculateElasticMatrix(VoigtSizeMatrixType& C, ConstitutiveLaw::Parameters& rValues) override;
157 
164  void CalculatePK2Stress(const Vector &rStrainVector, ConstitutiveLaw::StressVectorType &rStressVector, ConstitutiveLaw::Parameters &rValues) override;
165 
171  void CalculateCauchyGreenStrain(
173  ConstitutiveLaw::StrainVectorType& rStrainVector
174  ) override;
175 
180 
181 
186 
190  friend class Serializer;
191 
192  void save(Serializer& rSerializer) const override
193  {
195  }
196 
197  void load(Serializer& rSerializer) override
198  {
200  }
201 
202 
203 }; // Class AxisymElasticIsotropic
204 } // namespace Kratos.
Definition: axisym_elastic_isotropic.h:43
ProcessInfo ProcessInfoType
Definition: axisym_elastic_isotropic.h:49
SizeType GetStrainSize() const override
Definition: axisym_elastic_isotropic.h:95
std::size_t SizeType
Definition: axisym_elastic_isotropic.h:51
ConstitutiveLaw BaseType
Definition: axisym_elastic_isotropic.h:50
KRATOS_CLASS_POINTER_DEFINITION(AxisymElasticIsotropic)
Definition: constitutive_law.h:47
std::size_t SizeType
Definition: constitutive_law.h:82
Definition: elastic_isotropic_3d.h:53
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
#define KRATOS_SERIALIZE_SAVE_BASE_CLASS(Serializer, BaseType)
Definition: define.h:812
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
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