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_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 <string>
17 #include <iostream>
18 #include "includes/define.h"
19 
20 // Project includes
21 #include "includes/serializer.h"
23 
24 // Application includes
26 
27 namespace Kratos
28 {
29 /*
30  - structure of UMATs:
31  - Function to get stress, stiffness matrix
32  void umat(double* STRESS, double* STATEV, double** DDSDDE, double* SSE, double* SPD, double* SCD,
33  double* RPL, double* DDSDDT, double* DRPLDE, double* DRPLDT, double* STRAN, double* DSTRAN,
34  double* TIME, double* DTIME, double* TEMP, double* DTEMP, double* PREDEF, double* DPRED,
35  char* MATERL, int* NDI, int* NSHR, int* NTENS, int* NSTATV, double* PROPS,
36  int* NPROPS, double* COORDS, double** DROT, double* PNEWDT, double* CELENT, double** DFGRD0,
37  double** DFGRD1, int* NOEL, int* NPT, double* KSLAY, double* KSPT, int* KSTEP,
38  int* KINC);
39 
40 */
41 
42 typedef void(*pF_UMATMod) (double* STRESS, double* STATEV, double** DDSDDE, double* SSE, double* SPD, double* SCD,
43  double* rpl, double* ddsddt, double* drplde, double* drpldt, double* stran, double* dstran,
44  double* time, double* dtime, double* temp, double* dtemp, double* predef, double* dpred,
45  char* materl, int* ndi, int* nshr, int* ntens, int* nstatv, const double* props, int* nprops,
46  double* coords, double** drot, double* pnewdt, double* celent, double** dfgrd0,
47  double** dfgrd1, int* noel, int* npt, double* kslay, double* kspt, int* kstep,
48  int* kinc);
49 
50 
53 
56 
60 
64 
68 
72 
74 
77  class KRATOS_API(GEO_MECHANICS_APPLICATION) SmallStrainUMAT3DLaw: public ConstitutiveLaw
78  {
79  public:
80  // The process info type definition
82 
83  // The base class ConstitutiveLaw type definition
85 
87  using SizeType = std::size_t;
88 
90  static constexpr SizeType Dimension = N_DIM_3D;
91 
93  static constexpr SizeType VoigtSize = VOIGT_SIZE_3D;
94 
97 
98 
100  //@name Life Cycle
102 
103  SmallStrainUMAT3DLaw() = default;
104  ~SmallStrainUMAT3DLaw() override = default;
109 
113  ConstitutiveLaw::Pointer Clone() const override;
114 
119  void GetLawFeatures(Features& rFeatures) override;
120 
125  {
126  return Dimension;
127  }
128 
132  SizeType GetStrainSize() const override
133  {
134  return VoigtSize;
135  }
136 
142  {
143  return StrainMeasure_Infinitesimal;
144  }
145 
151  {
152  return StressMeasure_Cauchy;
153  }
154 
161  void CalculateMaterialResponsePK1(ConstitutiveLaw::Parameters & rValues) override;
162 
169  void CalculateMaterialResponsePK2(ConstitutiveLaw::Parameters & rValues) override;
170 
177  void CalculateMaterialResponseKirchhoff(ConstitutiveLaw::Parameters & rValues) override;
178 
185  void CalculateMaterialResponseCauchy(ConstitutiveLaw::Parameters & rValues) override;
186 
193  void FinalizeMaterialResponseCauchy (ConstitutiveLaw::Parameters & rValues) override;
194  void FinalizeMaterialResponsePK1 (ConstitutiveLaw::Parameters & rValues) override;
195  void FinalizeMaterialResponsePK2 (ConstitutiveLaw::Parameters & rValues) override;
196  void FinalizeMaterialResponseKirchhoff(ConstitutiveLaw::Parameters & rValues) override;
197 
205  double& CalculateValue(ConstitutiveLaw::Parameters& rParameterValues,
206  const Variable<double>& rThisVariable,
207  double& rValue) override;
208 
216  Vector& CalculateValue(ConstitutiveLaw::Parameters& rParameterValues,
217  const Variable<Vector>& rThisVariable,
218  Vector& rValue) override;
219 
227  Matrix& CalculateValue(ConstitutiveLaw::Parameters& rParameterValues,
228  const Variable<Matrix>& rThisVariable,
229  Matrix& rValue) override;
230 
232 
233 
234  // @brief This function provides the place to perform checks on the completeness of the input.
235  // @details It is designed to be called only once (or anyway, not often) typically at the beginning
236  // of the calculations, so to verify that nothing is missing from the input or that
237  // no common error is found.
238  int Check(const Properties& rMaterialProperties,
239  const GeometryType& rElementGeometry,
240  const ProcessInfo& rCurrentProcessInfo) const override;
241 
250  void InitializeMaterial(const Properties& rMaterialProperties,
251  const GeometryType& rElementGeometry,
252  const Vector& rShapeFunctionsValues) override;
253 
254 
259  void InitializeMaterialResponseCauchy (ConstitutiveLaw::Parameters& rValues) override;
260  void InitializeMaterialResponsePK1 (ConstitutiveLaw::Parameters& rValues) override;
261  void InitializeMaterialResponsePK2 (ConstitutiveLaw::Parameters& rValues) override;
262  void InitializeMaterialResponseKirchhoff(ConstitutiveLaw::Parameters& rValues) override;
263 
269  virtual void CalculateCauchyGreenStrain(ConstitutiveLaw::Parameters& rValues, Vector& rStrainVector);
270 
279  void ResetMaterial(const Properties& rMaterialProperties,
280  const GeometryType& rElementGeometry,
281  const Vector& rShapeFunctionsValues) override;
282 
284  double& GetValue( const Variable<double> &rThisVariable, double &rValue ) override;
285  int& GetValue( const Variable<int> &rThisVariable, int &rValue ) override;
286  Vector& GetValue( const Variable<Vector> &rThisVariable, Vector &rValue ) override;
287 
289  void SetValue( const Variable<double>& rVariable,
290  const double& rValue,
291  const ProcessInfo& rCurrentProcessInfo ) override;
292 
293  void SetValue( const Variable<Vector>& rVariable,
294  const Vector& rValue,
295  const ProcessInfo& rCurrentProcessInfo ) override;
296 
300 
304 
306  std::string Info() const override
307  {
308  return "SmallStrainUMAT3DLaw";
309  }
310 
312  void PrintInfo(std::ostream& rOStream) const override
313  {
314  rOStream << Info();
315  }
316 
318  void PrintData(std::ostream& rOStream) const override
319  {
320  rOStream << "SmallStrainUMAT3DLaw Data";
321  }
322 
326 
328 
329  protected:
334 
337 
338  double mMatrixD[VOIGT_SIZE_3D][VOIGT_SIZE_3D];
339 
343 
344 
348 
352  virtual void UpdateInternalDeltaStrainVector(ConstitutiveLaw::Parameters &rValues);
353  virtual void UpdateInternalStrainVectorFinalized(ConstitutiveLaw::Parameters &rValues);
354  virtual void SetExternalStressVector(Vector& rStressVector);
355  virtual void SetInternalStressVector(const Vector& rStressVector);
356  virtual void SetInternalStrainVector(const Vector& rStrainVector);
357  virtual void CopyConstitutiveMatrix(ConstitutiveLaw::Parameters &rValues, Matrix& rConstitutiveMatrix);
358 
359 
360  void CalculateConstitutiveMatrix(ConstitutiveLaw::Parameters &rValues, Matrix& rConstitutiveMatrix);
361  void CalculateStress(ConstitutiveLaw::Parameters &rValues, Vector& rStressVector);
362 
363  int GetStateVariableIndex(const Variable<double>& rThisVariable);
364 
368 
369 
373 
374 
376 
377  private:
380 
381 
385  pF_UMATMod pUserMod;
386 
387  bool mIsModelInitialized = false;
388  bool mIsUMATLoaded = false;
389 
390  std::vector<int> mProjectDirectory;
391 
392  Vector mStateVariables;
393  Vector mStateVariablesFinalized;
394 
395 
399 
400 
404 
405 
409 
410  // to load UMAT and functions
411  bool loadUMAT(const Properties& rMaterialProperties);
412  bool loadUMATWindows(const Properties& rMaterialProperties);
413  bool loadUMATLinux(const Properties& rMaterialProperties);
414 
415  // Set number of MaterialParameters
416  void CallUMAT( ConstitutiveLaw::Parameters &rValues);
417 
418  // Set state variables to the initial values
419  void ResetStateVariables(const Properties& rMaterialProperties);
420 
421 
425  friend class Serializer;
426 
427  void save(Serializer& rSerializer) const override
428  {
430  rSerializer.save("InitializedModel", mIsModelInitialized);
431  rSerializer.save("StressVectorFinalized", mStressVectorFinalized);
432  rSerializer.save("StrainVectorFinalized", mStrainVectorFinalized);
433  rSerializer.save("StateVariablesFinalized", mStateVariablesFinalized);
434  }
435 
436  void load(Serializer& rSerializer) override
437  {
439  rSerializer.load("InitializedModel", mIsModelInitialized);
440  rSerializer.load("StressVectorFinalized", mStressVectorFinalized);
441  rSerializer.load("StrainVectorFinalized", mStrainVectorFinalized);
442  rSerializer.load("StateVariablesFinalized", mStateVariablesFinalized);
443  }
444 
448 
449 
453 
455 
456  }; // Class SmallStrainUMAT3DLaw
457 
459 
462 
463 
467 
469 
471 
472 }
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
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
virtual bool & CalculateValue(Parameters &rParameterValues, const Variable< bool > &rThisVariable, bool &rValue)
Calculates the value of a specified variable (bool)
Definition: constitutive_law.cpp:370
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
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
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
Short class definition.
Definition: small_strain_umat_3D_law.hpp:78
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: small_strain_umat_3D_law.hpp:318
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: small_strain_umat_3D_law.hpp:312
SizeType GetStrainSize() const override
Voigt tensor size:
Definition: small_strain_umat_3D_law.hpp:132
~SmallStrainUMAT3DLaw() override=default
SmallStrainUMAT3DLaw & operator=(SmallStrainUMAT3DLaw &&)=delete
int GetStateVariableIndex(const Variable< double > &rThisVariable)
array_1d< double, VOIGT_SIZE_3D > mStressVectorFinalized
Definition: small_strain_umat_3D_law.hpp:333
StrainMeasure GetStrainMeasure() override
Returns the expected strain measure of this constitutive law (by default Green-Lagrange)
Definition: small_strain_umat_3D_law.hpp:141
SizeType WorkingSpaceDimension() override
Dimension of the law:
Definition: small_strain_umat_3D_law.hpp:124
array_1d< double, VOIGT_SIZE_3D > mStrainVectorFinalized
Definition: small_strain_umat_3D_law.hpp:336
SmallStrainUMAT3DLaw(SmallStrainUMAT3DLaw &&)=delete
KRATOS_CLASS_POINTER_DEFINITION(SmallStrainUMAT3DLaw)
Pointer definition of SmallStrainUMAT3DLaw.
array_1d< double, VOIGT_SIZE_3D > mDeltaStrainVector
Definition: small_strain_umat_3D_law.hpp:335
StressMeasure GetStressMeasure() override
Definition: small_strain_umat_3D_law.hpp:150
std::string Info() const override
Turn back information as a string.
Definition: small_strain_umat_3D_law.hpp:306
array_1d< double, VOIGT_SIZE_3D > mStressVector
Definition: small_strain_umat_3D_law.hpp:332
#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
void(* pF_UMATMod)(double *STRESS, double *STATEV, double **DDSDDE, double *SSE, double *SPD, double *SCD, double *rpl, double *ddsddt, double *drplde, double *drpldt, double *stran, double *dstran, double *time, double *dtime, double *temp, double *dtemp, double *predef, double *dpred, char *materl, int *ndi, int *nshr, int *ntens, int *nstatv, const double *props, int *nprops, double *coords, double **drot, double *pnewdt, double *celent, double **dfgrd0, double **dfgrd1, int *noel, int *npt, double *kslay, double *kspt, int *kstep, int *kinc)
Definition: small_strain_umat_3D_law.hpp:42
constexpr SizeType VOIGT_SIZE_3D
Definition: geo_mechanics_application_constants.h:41
def SetValue(entity, variable, value)
Definition: coupling_interface_data.py:256
time
Definition: face_heat.py:85
def load(f)
Definition: ode_solve.py:307
float temp
Definition: rotating_cone.py:85
Definition: constitutive_law.h:137
Definition: constitutive_law.h:189