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_3D_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) SmallStrainUMAT3DInterfaceLaw: public SmallStrainUMAT3DLaw
49  {
50  public:
51  // The process info type definition
53 
54  // The base class ConstitutiveLaw type definition
56 
58  using SizeType = std::size_t;
59 
61  static constexpr SizeType Dimension = N_DIM_3D;
62 
64  static constexpr SizeType VoigtSize = VOIGT_SIZE_3D_INTERFACE;
65 
68 
70  //@name Life Cycle
72 
76  ConstitutiveLaw::Pointer Clone() const override;
78  Vector& GetValue(const Variable<Vector> &rThisVariable, Vector &rValue) override;
79 
81  void SetValue(const Variable<Vector>& rVariable,
82  const Vector& rValue,
83  const ProcessInfo& rCurrentProcessInfo) override;
84 
89  {
90  return Dimension;
91  }
92 
96  SizeType GetStrainSize() const override
97  {
98  return VoigtSize;
99  }
100 
106  {
107  return StrainMeasure_Infinitesimal;
108  }
109 
115  {
116  return StressMeasure_Cauchy;
117  }
118 
119 
125  void CalculateCauchyGreenStrain(ConstitutiveLaw::Parameters& rValues, Vector& rStrainVector) override;
126 
130 
131 
135 
137  std::string Info() const override
138  {
139  return "SmallStrainUMAT3DInterfaceLaw";
140  }
141 
143  void PrintInfo(std::ostream& rOStream) const override
144  {
145  rOStream << Info();
146  }
147 
149  void PrintData(std::ostream& rOStream) const override
150  {
151  rOStream << "SmallStrainUMAT3DInterfaceLaw Data";
152  }
153 
154 
158 
159 
161 
162  protected:
165 
169 
173 
174 
178 
182  void UpdateInternalDeltaStrainVector(ConstitutiveLaw::Parameters &rValues) override;
183  void SetExternalStressVector(Vector& rStressVector) override;
184  void SetInternalStressVector(const Vector& rStressVector) override;
185  void SetInternalStrainVector(const Vector& rStrainVector) override;
186  void CopyConstitutiveMatrix(ConstitutiveLaw::Parameters &rValues, Matrix& rConstitutiveMatrix) override;
187 
188 
192 
193 
197 
198 
200 
201  private:
204 
205  indexStress3D getIndex3D(indexStress3DInterface index3D) const;
206 
210 
214 
215 
219 
220 
224 
225 
226 
230  friend class Serializer;
231 
232  void save(Serializer& rSerializer) const override
233  {
235  }
236 
237  void load(Serializer& rSerializer) override
238  {
240  }
241 
245 
246 
250 
252 
253  }; // Class SmallStrainUMAT3DLaw
254 
256 
259 
260 
264 
266 
268 
269 }
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_3D_interface_law.hpp:49
SizeType WorkingSpaceDimension() override
Dimension of the law:
Definition: small_strain_umat_3D_interface_law.hpp:88
StrainMeasure GetStrainMeasure() override
Returns the expected strain measure of this constitutive law (by default Green-Lagrange)
Definition: small_strain_umat_3D_interface_law.hpp:105
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: small_strain_umat_3D_interface_law.hpp:149
StressMeasure GetStressMeasure() override
Definition: small_strain_umat_3D_interface_law.hpp:114
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: small_strain_umat_3D_interface_law.hpp:143
KRATOS_CLASS_POINTER_DEFINITION(SmallStrainUMAT3DInterfaceLaw)
Pointer definition of SmallStrainUMAT3DInterfaceLaw.
SizeType GetStrainSize() const override
Voigt tensor size:
Definition: small_strain_umat_3D_interface_law.hpp:96
std::string Info() const override
Turn back information as a string.
Definition: small_strain_umat_3D_interface_law.hpp:137
Short class definition.
Definition: small_strain_umat_3D_law.hpp:78
virtual void SetValue(const Variable< bool > &rVariable, const bool &Value, const ProcessInfo &rCurrentProcessInfo)
Sets the value of a specified variable (boolean)
Definition: constitutive_law.cpp:279
virtual bool & GetValue(const Variable< bool > &rThisVariable, bool &rValue)
Returns the value of a specified variable (boolean)
Definition: constitutive_law.cpp:201
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_3D
Definition: geo_mechanics_application_constants.h:25
constexpr SizeType VOIGT_SIZE_3D_INTERFACE
Definition: geo_mechanics_application_constants.h:46
indexStress3D
Definition: geo_mechanics_application_constants.h:54
indexStress3DInterface
Definition: geo_mechanics_application_constants.h:83
def SetValue(entity, variable, value)
Definition: coupling_interface_data.py:256
def load(f)
Definition: ode_solve.py:307
Definition: constitutive_law.h:189