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.
thermal_elastic_isotropic_3d.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: Alejandro Cornejo
12 //
13 
14 # pragma once
15 
16 // System includes
17 
18 // External includes
19 
20 // Project includes
21 #include "custom_constitutive/elastic_isotropic_3d.h"
22 
23 namespace Kratos
24 {
27 
31 
35 
39 
43 
51 class KRATOS_API(CONSTITUTIVE_LAWS_APPLICATION) ThermalElasticIsotropic3D
52  : public ElasticIsotropic3D
53 {
54 public:
55 
58 
60 
63 
66 
71  {}
72 
76  ConstitutiveLaw::Pointer Clone() const override
77  {
78  return Kratos::make_shared<ThermalElasticIsotropic3D>(*this);
79  }
80 
81 
86 
91  : BaseType(rOther),
92  mReferenceTemperature(rOther.mReferenceTemperature)
93  {
94  }
95 
99 
103 
110  void CalculateMaterialResponsePK2(ConstitutiveLaw::Parameters &rValues) override;
111 
120  int Check(
121  const Properties& rMaterialProperties,
122  const GeometryType& rElementGeometry,
123  const ProcessInfo& rCurrentProcessInfo
124  ) const override;
125 
133  double& CalculateValue(
134  ConstitutiveLaw::Parameters& rParameterValues,
135  const Variable<double>& rThisVariable,
136  double& rValue
137  ) override;
138 
147  void InitializeMaterial(
148  const Properties &rMaterialProperties,
149  const GeometryType &rElementGeometry,
150  const Vector &rShapeFunctionsValues) override;
151 
152 
159  void CalculatePK2Stress(
160  const ConstitutiveLaw::StrainVectorType &rStrainVector,
161  ConstitutiveLaw::StressVectorType &rStressVector,
162  ConstitutiveLaw::Parameters &rValues) override;
166 
172  void CalculateElasticMatrix(
173  ConstitutiveLaw::VoigtSizeMatrixType& rConstitutiveMatrix,
175  ) override;
176 
184  virtual void SubstractThermalStrain(
185  ConstitutiveLaw::StrainVectorType &rStrainVector,
186  const double ReferenceTemperature,
187  ConstitutiveLaw::Parameters &rParameters,
188  const bool IsPlaneStrain = false);
189 
193 
197 
201 
203 
204 protected:
205 
208 
212 
216 
220 
222 
223 private:
224 
227 
231 
232  double mReferenceTemperature = 0.0;
233 
237 
242 
247 
252  double& GetReferenceTemperature()
253  {
254  return mReferenceTemperature;
255  }
256 
261  void SetReferenceTemperature(const double ToRefTemperature)
262  {
263  mReferenceTemperature = ToRefTemperature;
264  }
265 
269  friend class Serializer;
270 
271  void save(Serializer& rSerializer) const override
272  {
274  rSerializer.save("ReferenceTemperature", mReferenceTemperature);
275  }
276 
277  void load(Serializer& rSerializer) override
278  {
280  rSerializer.load("ReferenceTemperature", mReferenceTemperature);
281  }
282 
283 
284 }; // Class LinearPlaneStrain
285 } // namespace Kratos.
Definition: elastic_isotropic_3d.h:53
Geometry base class.
Definition: geometry.h:71
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
void load(std::string const &rTag, TDataType &rObject)
Definition: serializer.h:207
void save(std::string const &rTag, std::array< TDataType, TDataSize > const &rObject)
Definition: serializer.h:545
This class defines a Thermo dependant CL, including the addition of thermal expansion strains.
Definition: thermal_elastic_isotropic_3d.h:53
ThermalElasticIsotropic3D()
Default constructor.
Definition: thermal_elastic_isotropic_3d.h:70
ThermalElasticIsotropic3D(const ThermalElasticIsotropic3D &rOther)
Definition: thermal_elastic_isotropic_3d.h:90
KRATOS_CLASS_POINTER_DEFINITION(ThermalElasticIsotropic3D)
Counted pointer of LinearPlaneStrain.
ConstitutiveLaw::Pointer Clone() const override
Clone method.
Definition: thermal_elastic_isotropic_3d.h:76
~ThermalElasticIsotropic3D() override
Destructor.
Definition: thermal_elastic_isotropic_3d.h:85
#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
def load(f)
Definition: ode_solve.py:307
Definition: constitutive_law.h:189