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_udsm_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 the functions in PLAXIS UDSM:
31  - Function to get stress, stiffness matrix, attribute, number of state variables, ...
32  void User_Mod(int *IDTASK, int *IMOD, int *ISUNDR,
33  int *ISTEP, int *ITER, int *IEL, int *INT,
34  double *X, double *Y, double *Z,
35  double *TIME0, double *DTIME,
36  double *PROPS, double *SIG0, double *SWP0, double *STVAR0,
37  double *DEPS, double **D, double *BULKW,
38  double *SIG, double *SWP, double *STVAR, int *IPL,
39  int *NSTAT, int *NONSYM, int *ISTRSDEP, int *ITIMEDEP, int *ITANG,
40  int *IPRDIR, int *IPRJLEN, int *IABORT);
41 
42  void GetParamCount(int *IMOD, int *NPARAM);
43  void GetStateVarCount(int *IMOD, int *NSTVAR);
44 */
45 
46 typedef void(*pF_GetParamCount) (int *, int *);
47 typedef void(*pF_GetStateVarCount)(int *, int *);
48 typedef void(*pF_UserMod) (int *, int *, int *,
49  int *, int *, int *, int *,
50  double *, double *, double *,
51  double *, double *,
52  const double *, double *, double *, double *,
53  double *, double **, double *,
54  double *, double *, double *, int *,
55  int *, int *, int *, int *, int *,
56  int *, int *, int *);
57 
60 
63 
67 
71 
75 
79 
81 
84  class KRATOS_API(GEO_MECHANICS_APPLICATION) SmallStrainUDSM3DLaw: public ConstitutiveLaw
85  {
86  public:
87  // The base class ConstitutiveLaw type definition
89 
91  typedef std::size_t SizeType;
92 
94  static constexpr SizeType Dimension = N_DIM_3D;
95 
97  static constexpr SizeType VoigtSize = VOIGT_SIZE_3D;
98 
101 
103  //@name Life Cycle
105 
106  //----------------------------------------------------------------------------------------
110  SmallStrainUDSM3DLaw() = default;
111 
115  ConstitutiveLaw::Pointer Clone() const override;
116 
121 
125  ~SmallStrainUDSM3DLaw() override = default;
126 
127  // Assignment operator:
129 
134  void GetLawFeatures(Features& rFeatures) override;
135 
140  {
141  return Dimension;
142  }
143 
147  SizeType GetStrainSize() const override
148  {
149  return VoigtSize;
150  }
151 
157  {
158  return StrainMeasure_Infinitesimal;
159  }
160 
166  {
167  return StressMeasure_Cauchy;
168  }
169 
176  void CalculateMaterialResponsePK1(ConstitutiveLaw::Parameters & rValues) override;
177 
184  void CalculateMaterialResponsePK2(ConstitutiveLaw::Parameters & rValues) override;
185 
192  void CalculateMaterialResponseKirchhoff(ConstitutiveLaw::Parameters & rValues) override;
193 
200  void CalculateMaterialResponseCauchy(ConstitutiveLaw::Parameters & rValues) override;
201 
208  void FinalizeMaterialResponseCauchy (ConstitutiveLaw::Parameters & rValues) override;
209  void FinalizeMaterialResponsePK1 (ConstitutiveLaw::Parameters & rValues) override;
210  void FinalizeMaterialResponsePK2 (ConstitutiveLaw::Parameters & rValues) override;
211  void FinalizeMaterialResponseKirchhoff(ConstitutiveLaw::Parameters & rValues) override;
212 
220  double& CalculateValue(ConstitutiveLaw::Parameters& rParameterValues,
221  const Variable<double>& rThisVariable,
222  double& rValue) override;
223 
231  Vector& CalculateValue(ConstitutiveLaw::Parameters& rParameterValues,
232  const Variable<Vector>& rThisVariable,
233  Vector& rValue) override;
234 
242  Matrix& CalculateValue(ConstitutiveLaw::Parameters& rParameterValues,
243  const Variable<Matrix>& rThisVariable,
244  Matrix& rValue) override;
245 
247 
248  // @brief This function provides the place to perform checks on the completeness of the input.
249  // @details It is designed to be called only once (or anyway, not often) typically at the beginning
250  // of the calculations, so to verify that nothing is missing from the input or that
251  // no common error is found.
252  int Check(const Properties& rMaterialProperties,
253  const GeometryType& rElementGeometry,
254  const ProcessInfo& rCurrentProcessInfo) const override;
255 
264  void InitializeMaterial(const Properties& rMaterialProperties,
265  const GeometryType& rElementGeometry,
266  const Vector& rShapeFunctionsValues) override;
267 
268 
273  void InitializeMaterialResponseCauchy (ConstitutiveLaw::Parameters& rValues) override;
274  void InitializeMaterialResponsePK1 (ConstitutiveLaw::Parameters& rValues) override;
275  void InitializeMaterialResponsePK2 (ConstitutiveLaw::Parameters& rValues) override;
276  void InitializeMaterialResponseKirchhoff(ConstitutiveLaw::Parameters& rValues) override;
277 
283  virtual void CalculateCauchyGreenStrain(ConstitutiveLaw::Parameters& rValues, Vector& rStrainVector);
284 
293  void ResetMaterial(const Properties& rMaterialProperties,
294  const GeometryType& rElementGeometry,
295  const Vector& rShapeFunctionsValues) override;
296 
297  double& GetValue( const Variable<double> &rThisVariable, double &rValue ) override;
298  int& GetValue( const Variable<int> &rThisVariable, int &rValue ) override;
299  Vector& GetValue( const Variable<Vector> &rThisVariable, Vector &rValue ) override;
300 
301  void SetValue(const Variable<double>& rVariable,
302  const double& rValue,
303  const ProcessInfo& rCurrentProcessInfo) override;
304 
305  void SetValue(const Variable<Vector>& rVariable,
306  const Vector& rValue,
307  const ProcessInfo& rCurrentProcessInfo) override;
311 
315 
317  std::string Info() const override
318  {
319  return "SmallStrainUDSM3DLaw";
320  }
321 
323  void PrintInfo(std::ostream& rOStream) const override
324  {
325  rOStream << Info();
326  }
327 
329  void PrintData(std::ostream& rOStream) const override
330  {
331  rOStream << "SmallStrainUDSM3DLaw Data";
332  }
333 
337 
339 
340  protected:
343  enum IDTASK: int{
344  INITIALISATION = 1,
349  MATRIX_ELASTIC
350  };
351 
352  enum ATTRIBUTE: int {
356  USE_TANGENT_MATRIX
357  };
358 
364 
367 
368  double mMatrixD[VOIGT_SIZE_3D][VOIGT_SIZE_3D];
369 
373 
377 
381  virtual void UpdateInternalDeltaStrainVector(ConstitutiveLaw::Parameters &rValues);
382  virtual void UpdateInternalStrainVectorFinalized(ConstitutiveLaw::Parameters &rValues);
383  virtual void SetExternalStressVector(Vector& rStressVector);
384  virtual void SetInternalStressVector(const Vector& rStressVector);
385  virtual void SetInternalStrainVector(const Vector& rStrainVector);
386  virtual void CopyConstitutiveMatrix(ConstitutiveLaw::Parameters &rValues, Matrix& rConstitutiveMatrix);
387 
388  void CalculateConstitutiveMatrix(ConstitutiveLaw::Parameters &rValues, Matrix& rConstitutiveMatrix);
389  void CalculateStress(ConstitutiveLaw::Parameters &rValues, Vector& rStressVector);
390 
391  // returns 1 if the stiffness matrix of the material is non-symmetric
392  int getIsNonSymmetric() {return mAttributes[IS_NON_SYMMETRIC];}
393 
394  // returns 1 if the stiffness matrix of the material is stress dependent
395  int getIsStressDependent() {return mAttributes[IS_STRESS_DEPENDENT];}
396 
397  // returns 1 if material is time dependent
398  int getIsTimeDependent() {return mAttributes[IS_TIME_DEPENDENT];}
399 
400  // returns 1 if the stiffness matrix of the material is tangential
401  int getUseTangentMatrix() {return mAttributes[USE_TANGENT_MATRIX];}
402 
406 
410 
412 
413  private:
416 
420  pF_GetParamCount pGetParamCount;
421  pF_GetStateVarCount pGetStateVarCount;
422  pF_UserMod pUserMod;
423 
424  bool mIsModelInitialized = false;
425  bool mIsUDSMLoaded = false;
426 
427  //AttributesUDSM mAttributes;
428  array_1d<int, 4> mAttributes;
429 
430  std::vector<int> mProjectDirectory;
431 
432  Vector mStateVariables;
433  Vector mStateVariablesFinalized;
434 
438 
442 
446 
447  // to load UDSM and functions
448  bool loadUDSM(const Properties& rMaterialProperties);
449  bool loadUDSMWindows(const Properties& rMaterialProperties);
450  bool loadUDSMLinux(const Properties& rMaterialProperties);
451 
452  // Set number of MaterialParameters
453  void CallUDSM(int *IDTask, ConstitutiveLaw::Parameters &rValues);
454 
455  // Set state variables to zero
456  void ResetStateVariables(const Properties& rMaterialProperties);
457 
458  // set number of StateVariables
459  void SetAttributes(const Properties& rMaterialProperties);
460 
461  // set number of StateVariables
462  int GetNumberOfStateVariablesFromUDSM(const Properties& rMaterialProperties);
463 
464  // get number of MaterialParameters
465  SizeType GetNumberOfMaterialParametersFromUDSM(const Properties& rMaterialProperties);
466  int GetStateVariableIndex(const Variable<double>& rThisVariable) const;
467 
471  friend class Serializer;
472 
473  void save(Serializer& rSerializer) const override
474  {
476  rSerializer.save("InitializedModel", mIsModelInitialized);
477  rSerializer.save("Attributes", mAttributes);
478  rSerializer.save("StressVectorFinalized", mStressVectorFinalized);
479  rSerializer.save("StrainVectorFinalized", mStrainVectorFinalized);
480  rSerializer.save("StateVariablesFinalized", mStateVariablesFinalized);
481  }
482 
483  void load(Serializer& rSerializer) override
484  {
486  rSerializer.load("InitializedModel", mIsModelInitialized);
487  rSerializer.load("Attributes", mAttributes);
488  rSerializer.load("StressVectorFinalized", mStressVectorFinalized);
489  rSerializer.load("StrainVectorFinalized", mStrainVectorFinalized);
490  rSerializer.load("StateVariablesFinalized", mStateVariablesFinalized);
491  }
492 
496 
500 
502 
503  }; // Class SmallStrainUDSM3DLaw
504 
506 
509 
513 
515 
517 
518 }
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
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_udsm_3D_law.hpp:85
StrainMeasure GetStrainMeasure() override
Returns the expected strain measure of this constitutive law (by default Green-Lagrange)
Definition: small_strain_udsm_3D_law.hpp:156
int getIsStressDependent()
Definition: small_strain_udsm_3D_law.hpp:395
int getIsTimeDependent()
Definition: small_strain_udsm_3D_law.hpp:398
~SmallStrainUDSM3DLaw() override=default
Destructor.
ATTRIBUTE
Definition: small_strain_udsm_3D_law.hpp:352
@ IS_TIME_DEPENDENT
Definition: small_strain_udsm_3D_law.hpp:355
@ IS_STRESS_DEPENDENT
Definition: small_strain_udsm_3D_law.hpp:354
@ IS_NON_SYMMETRIC
Definition: small_strain_udsm_3D_law.hpp:353
IDTASK
Definition: small_strain_udsm_3D_law.hpp:343
@ MATRIX_ELASTO_PLASTIC
Definition: small_strain_udsm_3D_law.hpp:346
@ STRESS_CALCULATION
Definition: small_strain_udsm_3D_law.hpp:345
@ ATTRIBUTES
Definition: small_strain_udsm_3D_law.hpp:348
@ NUMBER_OF_STATE_VARIABLES
Definition: small_strain_udsm_3D_law.hpp:347
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: small_strain_udsm_3D_law.hpp:323
int getUseTangentMatrix()
Definition: small_strain_udsm_3D_law.hpp:401
SmallStrainUDSM3DLaw()=default
Default constructor.
array_1d< double, VOIGT_SIZE_3D > mStressVectorFinalized
Definition: small_strain_udsm_3D_law.hpp:363
SizeType WorkingSpaceDimension() override
Dimension of the law:
Definition: small_strain_udsm_3D_law.hpp:139
array_1d< double, VOIGT_SIZE_3D > mDeltaStrainVector
Definition: small_strain_udsm_3D_law.hpp:365
std::size_t SizeType
The size type definition.
Definition: small_strain_udsm_3D_law.hpp:91
std::string Info() const override
Turn back information as a string.
Definition: small_strain_udsm_3D_law.hpp:317
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: small_strain_udsm_3D_law.hpp:329
StressMeasure GetStressMeasure() override
Definition: small_strain_udsm_3D_law.hpp:165
int getIsNonSymmetric()
Definition: small_strain_udsm_3D_law.hpp:392
array_1d< double, VOIGT_SIZE_3D > mStrainVectorFinalized
Definition: small_strain_udsm_3D_law.hpp:366
KRATOS_CLASS_POINTER_DEFINITION(SmallStrainUDSM3DLaw)
Pointer definition of SmallStrainUDSM3DLaw.
ConstitutiveLaw BaseType
Definition: small_strain_udsm_3D_law.hpp:88
array_1d< double, VOIGT_SIZE_3D > mStressVector
Definition: small_strain_udsm_3D_law.hpp:362
SizeType GetStrainSize() const override
Voigt tensor size:
Definition: small_strain_udsm_3D_law.hpp:147
#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
void(* pF_GetStateVarCount)(int *, int *)
Definition: small_strain_udsm_3D_law.hpp:47
constexpr SizeType N_DIM_3D
Definition: geo_mechanics_application_constants.h:25
std::size_t SizeType
The definition of the size type.
Definition: mortar_classes.h:43
void(* pF_UserMod)(int *, int *, int *, int *, int *, int *, int *, double *, double *, double *, double *, double *, const double *, double *, double *, double *, double *, double **, double *, double *, double *, double *, int *, int *, int *, int *, int *, int *, int *, int *, int *)
Definition: small_strain_udsm_3D_law.hpp:48
constexpr SizeType VOIGT_SIZE_3D
Definition: geo_mechanics_application_constants.h:41
void(* pF_GetParamCount)(int *, int *)
Definition: small_strain_udsm_3D_law.hpp:46
def SetValue(entity, variable, value)
Definition: coupling_interface_data.py:256
def load(f)
Definition: ode_solve.py:307
Definition: constitutive_law.h:137
Definition: constitutive_law.h:189