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.
user_provided_linear_elastic_law.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: Ruben Zorrilla
10 // Riccardo Rossi
11 //
12 
13 #pragma once
14 
15 // System includes
16 
17 // External includes
18 
19 // Project includes
21 
22 namespace Kratos
23 {
24 
27 
31 
35 
39 
43 
53 template <unsigned int TDim>
54 class KRATOS_API(STRUCTURAL_MECHANICS_APPLICATION) UserProvidedLinearElasticLaw
55  : public ConstitutiveLaw
56 {
57 public:
58 
61 
64 
67 
69  typedef std::size_t SizeType;
70 
72  static constexpr SizeType Dimension = TDim;
73 
75  static constexpr SizeType StrainSize = (TDim * 3) - 3;
76 
79 
83 
88 
92  ConstitutiveLaw::Pointer Clone() const override;
93 
98 
103 
107 
111 
116  void GetLawFeatures(Features& rFeatures) override;
117 
122  {
123  return Dimension;
124  };
125 
129  SizeType GetStrainSize() const override
130  {
131  return StrainSize;
132  };
133 
139  {
140  return StrainMeasure_Infinitesimal;
141  }
142 
148  {
149  return StressMeasure_Cauchy;
150  }
151 
158  void CalculateMaterialResponsePK1 (ConstitutiveLaw::Parameters & rValues) override;
159 
166  void CalculateMaterialResponsePK2 (ConstitutiveLaw::Parameters & rValues) override;
167 
174  void CalculateMaterialResponseKirchhoff (ConstitutiveLaw::Parameters & rValues) override;
175 
182  void CalculateMaterialResponseCauchy (ConstitutiveLaw::Parameters & rValues) override;
183 
189  void InitializeMaterialResponsePK1 (Parameters& rValues) override {};
190 
196  void InitializeMaterialResponsePK2 (Parameters& rValues) override {};
197 
203  void InitializeMaterialResponseKirchhoff (Parameters& rValues) override {};
204 
210  void InitializeMaterialResponseCauchy (Parameters& rValues) override {};
211 
217  void FinalizeMaterialResponsePK1 (Parameters& rValues) override {};
218 
224  void FinalizeMaterialResponsePK2 (Parameters& rValues) override {};
225 
231  void FinalizeMaterialResponseKirchhoff (Parameters& rValues) override {};
232 
239  void FinalizeMaterialResponseCauchy (Parameters& rValues) override {};
240 
248  double& CalculateValue(
249  ConstitutiveLaw::Parameters& rParameterValues,
250  const Variable<double>& rThisVariable,
251  double& rValue) override;
252 
260  Vector& CalculateValue(
261  ConstitutiveLaw::Parameters& rParameterValues,
262  const Variable<Vector>& rThisVariable,
263  Vector& rValue) override;
264 
265  // /**
266  // * @brief It calculates the value of a specified variable (Vector case)
267  // * @param rParameterValues the needed parameters for the CL calculation
268  // * @param rThisVariable the variable to be returned
269  // * @param rValue a reference to the returned value
270  // * @return rValue output: the value of the specified variable
271  // */
272  // ConstitutiveLaw::StrainVectorType& CalculateValue(
273  // ConstitutiveLaw::Parameters& rParameterValues,
274  // const Variable<StrainVectorType>& rThisVariable,
275  // ConstitutiveLaw::StrainVectorType& rValue) override;
276 
277  // /**
278  // * @brief It calculates the value of a specified variable (Vector case)
279  // * @param rParameterValues the needed parameters for the CL calculation
280  // * @param rThisVariable the variable to be returned
281  // * @param rValue a reference to the returned value
282  // * @return rValue output: the value of the specified variable
283  // */
284  // ConstitutiveLaw::StressVectorType& CalculateValue(
285  // ConstitutiveLaw::Parameters& rParameterValues,
286  // const Variable<StressVectorType>& rThisVariable,
287  // ConstitutiveLaw::StressVectorType& rValue) override;
288 
296  Matrix& CalculateValue(
297  ConstitutiveLaw::Parameters& rParameterValues,
298  const Variable<Matrix>& rThisVariable,
299  Matrix& rValue) override;
300 
301  // /**
302  // * @brief It calculates the value of a specified variable (Matrix case)
303  // * @param rParameterValues the needed parameters for the CL calculation
304  // * @param rThisVariable the variable to be returned
305  // * @param rValue a reference to the returned value
306  // * @return rValue output: the value of the specified variable
307  // */
308  // ConstitutiveLaw::VoigtSizeMatrixType& CalculateValue(
309  // ConstitutiveLaw::Parameters& rParameterValues,
310  // const Variable<VoigtSizeMatrixType>& rThisVariable,
311  // ConstitutiveLaw::VoigtSizeMatrixType& rValue) override;
312 
313  // /**
314  // * @brief It calculates the value of a specified variable (Matrix case)
315  // * @param rParameterValues the needed parameters for the CL calculation
316  // * @param rThisVariable the variable to be returned
317  // * @param rValue a reference to the returned value
318  // * @return rValue output: the value of the specified variable
319  // */
320  // ConstitutiveLaw::DeformationGradientMatrixType& CalculateValue(
321  // ConstitutiveLaw::Parameters& rParameterValues,
322  // const Variable<DeformationGradientMatrixType>& rThisVariable,
323  // ConstitutiveLaw::DeformationGradientMatrixType& rValue) override;
324 
333  int Check(
334  const Properties& rMaterialProperties,
335  const GeometryType& rElementGeometry,
336  const ProcessInfo& rCurrentProcessInfo) const override;
337 
338 protected:
339 
342 
346 
350 
354 
360  virtual void CalculateElasticMatrix(
361  ConstitutiveLaw::VoigtSizeMatrixType& rConstitutiveMatrix,
363  );
364 
371  virtual void CalculatePK2Stress(
372  const ConstitutiveLaw::StrainVectorType& rStrainVector,
373  ConstitutiveLaw::StressVectorType& rStressVector,
375  );
376 
382  virtual void CalculateGreenLagrangeStrainVector(
384  ConstitutiveLaw::StrainVectorType& rStrainVector
385  );
386 
388 
389 private:
390 
393 
397 
401 
405 
410 
414  friend class Serializer;
415 
416  void save(Serializer& rSerializer) const override
417  {
419  }
420 
421  void load(Serializer& rSerializer) override
422  {
424  }
425 
426 
427 }; // Class UserProvidedLinearElasticLaw
428 } // namespace Kratos.
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
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
This class defines a linear elastic law with user provided constitutive tensor.
Definition: user_provided_linear_elastic_law.h:56
void FinalizeMaterialResponsePK2(Parameters &rValues) override
Definition: user_provided_linear_elastic_law.h:224
void FinalizeMaterialResponsePK1(Parameters &rValues) override
Definition: user_provided_linear_elastic_law.h:217
void FinalizeMaterialResponseKirchhoff(Parameters &rValues) override
Definition: user_provided_linear_elastic_law.h:231
SizeType GetStrainSize() const override
Voigt tensor size:
Definition: user_provided_linear_elastic_law.h:129
StressMeasure GetStressMeasure() override
Returns the stress measure of this constitutive law (by default 2st Piola-Kirchhoff stress in voigt n...
Definition: user_provided_linear_elastic_law.h:147
SizeType WorkingSpaceDimension() override
Dimension of the law:
Definition: user_provided_linear_elastic_law.h:121
void InitializeMaterialResponseKirchhoff(Parameters &rValues) override
Definition: user_provided_linear_elastic_law.h:203
ConstitutiveLaw BaseType
The base class ConstitutiveLaw type definition.
Definition: user_provided_linear_elastic_law.h:66
StrainMeasure GetStrainMeasure() override
Returns the expected strain measure of this constitutive law (by default Green-Lagrange)
Definition: user_provided_linear_elastic_law.h:138
KRATOS_CLASS_POINTER_DEFINITION(UserProvidedLinearElasticLaw)
Counted pointer of UserProvidedLinearElasticLaw.
void FinalizeMaterialResponseCauchy(Parameters &rValues) override
Definition: user_provided_linear_elastic_law.h:239
void InitializeMaterialResponsePK2(Parameters &rValues) override
Definition: user_provided_linear_elastic_law.h:196
void InitializeMaterialResponseCauchy(Parameters &rValues) override
Definition: user_provided_linear_elastic_law.h:210
ProcessInfo ProcessInfoType
The process info type definition.
Definition: user_provided_linear_elastic_law.h:63
std::size_t SizeType
The size type definition.
Definition: user_provided_linear_elastic_law.h:69
void InitializeMaterialResponsePK1(Parameters &rValues) override
Definition: user_provided_linear_elastic_law.h:189
#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:137
Definition: constitutive_law.h:189