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_orthotropic_2D_law.h
Go to the documentation of this file.
1 // KRATOS ___ _ _ _ _ _ __ _
2 // / __\___ _ __ ___| |_(_) |_ _ _| |_(_)_ _____ / / __ ___ _____ /_\ _ __ _ __
3 // / / / _ \| '_ \/ __| __| | __| | | | __| \ \ / / _ \/ / / _` \ \ /\ / / __| //_\\| '_ \| '_ |
4 // / /__| (_) | | | \__ \ |_| | |_| |_| | |_| |\ V / __/ /__| (_| |\ V V /\__ \/ _ \ |_) | |_) |
5 // \____/\___/|_| |_|___/\__|_|\__|\__,_|\__|_| \_/ \___\____/\__,_| \_/\_/ |___/\_/ \_/ .__/| .__/
6 // |_| |_|
7 //
8 // License: BSD License
9 // license: structural_mechanics_application/license.txt
10 //
11 // Main authors: Peter Wilson
12 // Contact: A.Winterstein [at] tum.de
13 //
14 
15 #pragma once
16 
17 // System includes
18 
19 // External includes
20 
21 // Project includes
23 
24 namespace Kratos
25 {
34  class KRATOS_API(CONSTITUTIVE_LAWS_APPLICATION) LinearElasticOrthotropic2DLaw : public ConstitutiveLaw
35  {
36  public:
42  typedef std::size_t SizeType;
48 
57 
62  ConstitutiveLaw::Pointer Clone() const override;
63 
68 
73 
84  SizeType GetStrainSize() const override
85  {
86  return 3;
87  };
88 
93  void GetLawFeatures(Features& rFeatures) override;
94 
101  void CalculateMaterialResponsePK2(ConstitutiveLaw::Parameters & rValues) override;
102 
109  bool& GetValue(const Variable<bool>& rThisVariable, bool& rValue) override;
110 
120  int Check(
121  const Properties& rMaterialProperties,
122  const GeometryType& rElementGeometry,
123  const ProcessInfo& rCurrentProcessInfo) const override;
124 
125  protected:
126 
139 
144  void CalculateGreenLagrangeStrain(const Matrix & rRightCauchyGreen,
145  Vector& rStrainVector);
146 
153  virtual void CalculateStress(const Vector &rStrainVector,
154  const Matrix &rConstitutiveMatrix,
155  Vector& rStressVector);
156 
165  virtual void CalculateLinearElasticMatrix(Matrix& rConstitutiveMatrix,
166  const Properties& rMaterialProperties);
167 
174  bool CheckParameters(ConstitutiveLaw::Parameters& rValues);
175 
176  private:
177 
190 
195 
199  friend class Serializer;
200 
201  void save(Serializer& rSerializer) const override
202  {
204  }
205 
206  void load(Serializer& rSerializer) override
207  {
209  }
210 
212  }; // Class LinearElasticOrthotropic2DLaw
213 } // namespace Kratos.
Definition: constitutive_law.h:47
std::size_t SizeType
Definition: constitutive_law.h:82
Geometry base class.
Definition: geometry.h:71
Definition: linear_elastic_orthotropic_2D_law.h:35
ConstitutiveLaw BaseType
Definition: linear_elastic_orthotropic_2D_law.h:41
KRATOS_CLASS_POINTER_DEFINITION(LinearElasticOrthotropic2DLaw)
SizeType GetStrainSize() const override
Definition: linear_elastic_orthotropic_2D_law.h:84
std::size_t SizeType
Definition: linear_elastic_orthotropic_2D_law.h:42
ProcessInfo ProcessInfoType
Definition: linear_elastic_orthotropic_2D_law.h:40
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
Properties encapsulates data shared by different Elements or Conditions. It can store any type of dat...
Definition: properties.h:69
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
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
def load(f)
Definition: ode_solve.py:307
Definition: constitutive_law.h:189