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.
small_strain_umat_2D_interface_law.hpp
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 // System includes
16 #include "includes/define.h"
17 
18 // Project includes
20 
21 namespace Kratos
22 {
25 
28 
32 
36 
40 
44 
46 
48  class KRATOS_API(GEO_MECHANICS_APPLICATION) SmallStrainUMAT2DInterfaceLaw: public SmallStrainUMAT3DLaw
49  {
50  public:
51  // The base class ConstitutiveLaw type definition
53 
55  using SizeType = std::size_t;
56 
58  static constexpr SizeType Dimension = N_DIM_2D;
59 
61  static constexpr SizeType VoigtSize = VOIGT_SIZE_2D_INTERFACE;
62 
65 
66 
68  //@name Life Cycle
70 
74  ConstitutiveLaw::Pointer Clone() const override;
75 
76  Vector& GetValue( const Variable<Vector> &rThisVariable, Vector &rValue ) override;
77 
78  void SetValue(const Variable<Vector>& rVariable,
79  const Vector& rValue,
80  const ProcessInfo& rCurrentProcessInfo ) override;
81 
86  {
87  return Dimension;
88  }
89 
93  SizeType GetStrainSize() const override
94  {
95  return VoigtSize;
96  }
97 
103  {
104  return StrainMeasure_Infinitesimal;
105  }
106 
112  {
113  return StressMeasure_Cauchy;
114  }
115 
116 
122  void CalculateCauchyGreenStrain(ConstitutiveLaw::Parameters& rValues, Vector& rStrainVector) override;
123 
127 
128 
132 
134  std::string Info() const override
135  {
136  return "SmallStrainUMAT2DInterfaceLaw";
137  }
138 
140  void PrintInfo(std::ostream& rOStream) const override
141  {
142  rOStream << Info();
143  }
144 
146  void PrintData(std::ostream& rOStream) const override
147  {
148  rOStream << "SmallStrainUMAT2DInterfaceLaw Data";
149  }
150 
154 
156 
157  protected:
160 
164 
168 
169 
173 
177  void UpdateInternalDeltaStrainVector(ConstitutiveLaw::Parameters &rValues) override;
178  void SetExternalStressVector(Vector& rStressVector) override;
179  void SetInternalStressVector(const Vector& rStressVector) override;
180  void SetInternalStrainVector(const Vector& rStrainVector) override;
181  void CopyConstitutiveMatrix(ConstitutiveLaw::Parameters &rValues, Matrix& rConstitutiveMatrix) override;
182 
186 
190 
192 
193  private:
196 
197  indexStress3D getIndex3D(indexStress2DInterface index2D) const;
198 
202 
206 
207 
211 
212 
216 
220  friend class Serializer;
221 
222  void save(Serializer& rSerializer) const override
223  {
225  }
226 
227  void load(Serializer& rSerializer) override
228  {
230  }
231 
235 
239 
241 
242  }; // Class SmallStrainUMAT3DLaw
243 
245 
248 
252 
254 
256 
257 }
std::string Info() const override
Turn back information as a string.
Definition: periodic_interface_process.hpp:93
Definition: constitutive_law.h:47
StrainMeasure
Definition: constitutive_law.h:52
StressMeasure
Definition: constitutive_law.h:69
std::size_t SizeType
Definition: constitutive_law.h:82
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
Short class definition.
Definition: small_strain_umat_2D_interface_law.hpp:49
std::string Info() const override
Turn back information as a string.
Definition: small_strain_umat_2D_interface_law.hpp:134
StrainMeasure GetStrainMeasure() override
Returns the expected strain measure of this constitutive law (by default Green-Lagrange)
Definition: small_strain_umat_2D_interface_law.hpp:102
KRATOS_CLASS_POINTER_DEFINITION(SmallStrainUMAT2DInterfaceLaw)
Pointer definition of SmallStrainUMAT2DInterfaceLaw.
SizeType GetStrainSize() const override
Voigt tensor size:
Definition: small_strain_umat_2D_interface_law.hpp:93
SizeType WorkingSpaceDimension() override
Dimension of the law:
Definition: small_strain_umat_2D_interface_law.hpp:85
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: small_strain_umat_2D_interface_law.hpp:146
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: small_strain_umat_2D_interface_law.hpp:140
StressMeasure GetStressMeasure() override
Definition: small_strain_umat_2D_interface_law.hpp:111
Short class definition.
Definition: small_strain_umat_3D_law.hpp:78
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
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
constexpr SizeType N_DIM_2D
Definition: geo_mechanics_application_constants.h:26
constexpr SizeType VOIGT_SIZE_2D_INTERFACE
Definition: geo_mechanics_application_constants.h:45
indexStress2DInterface
Definition: geo_mechanics_application_constants.h:79
indexStress3D
Definition: geo_mechanics_application_constants.h:54
def SetValue(entity, variable, value)
Definition: coupling_interface_data.py:256
def load(f)
Definition: ode_solve.py:307
Definition: constitutive_law.h:189