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_model.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosUmatApplication $
3 // Created by: $Author: LMonforte $
4 // Last modified by: $Co-Author: $
5 // Date: $Date: October 2017 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_SMALL_STRAIN_UMAT_MODEL_H_INCLUDED )
11 #define KRATOS_SMALL_STRAIN_UMAT_MODEL_H_INCLUDED
12 
13 // System includes
14 #include <string>
15 #include <iostream>
16 
17 // External includes
18 
19 // Project includes
22 
23 namespace Kratos
24 {
27 
30 
34 
38 
42 
46 
48 
50  class KRATOS_API(UMAT_APPLICATION) SmallStrainUmatModel : public ConstitutiveModel
51  {
52  protected:
53 
55  {
56  private:
57 
58  Flags* mpState;
59  const ModelDataType* mpModelData;
60 
61  public:
62 
65 
66  //Set Data Pointers
67  void SetState (Flags& rState) {mpState = &rState;};
68  void SetModelData (const ModelDataType& rModelData) {mpModelData = &rModelData;};
69 
70  //Get Data Pointers
71  const ModelDataType& GetModelData () const {return *mpModelData;};
72  const MaterialDataType& GetMaterialParameters () const {return mpModelData->GetMaterialParameters();};
73 
74  //Get non const Data
75  Flags& State () {return *mpState;};
76 
77  //Get const Data
78  const Flags& GetState () const {return *mpState;};
79 
80  };
81 
82 
83  public:
84 
88 
89 
92 
96 
99 
102 
104  virtual ConstitutiveModel::Pointer Clone() const override;
105 
108 
110  virtual ~SmallStrainUmatModel();
111 
112 
116 
117 
121 
122 
126  virtual void InitializeModel(ModelDataType& rValues) override;
127 
128 
132  virtual void FinalizeModel(ModelDataType& rValues) override;
133 
134 
138  virtual void CalculateStrainEnergy(ModelDataType& rValues, double& rDensityFunction) override;
139 
140 
144  virtual void CalculateStressTensor(ModelDataType& rValues, MatrixType& rStressMatrix) override;
145 
146 
147 
148 
152  virtual void CalculateConstitutiveTensor(ModelDataType& rValues, Matrix& rConstitutiveMatrix) override;
153 
154 
158  virtual void CalculateStressAndConstitutiveTensors(ModelDataType& rValues, MatrixType& rStressMatrix, Matrix& rConstitutiveMatrix) override;
159 
160 
164  virtual int Check(const Properties& rMaterialProperties, const ProcessInfo& rCurrentProcessInfo) override;
165 
169 
170  virtual void SetValue(const Variable<Vector>& rThisVariable, const Vector& rValue,
171  const ProcessInfo& rCurrentProcessInfo ) override
172  {
173  KRATOS_TRY
174 
175  // A method to compute the initial linear strain from the stress is needed
176  //if(rThisVariable == INITIAL_STRESS_VECTOR)
177 
178  // A method to compute the initial linear strain from the stress is needed
179  // if(rThisVariable == INITIAL_STRAIN_VECTOR){
180  // this->mHistoryVector = rValue;
181  // }
182 
183  KRATOS_CATCH(" ")
184  }
185 
186 
187  virtual void SetValue(const Variable<Matrix>& rThisVariable, const Matrix& rValue,
188  const ProcessInfo& rCurrentProcessInfo ) override
189  {
190  KRATOS_TRY
191 
192  // A method to compute the initial linear strain from the stress is needed
193  //if(rThisVariable == INITIAL_STRESS_VECTOR)
194 
195  // A method to compute the initial linear strain from the stress is needed
196  // if(rThisVariable == INITIAL_STRAIN_VECTOR){
197  // this->mHistoryVector = rValue;
198  // }
199 
200  KRATOS_CATCH(" ")
201  }
202 
208  virtual void GetDomainVariablesList(std::vector<Variable<double> >& rScalarVariables,
209  std::vector<Variable<array_1d<double,3> > >& rComponentVariables) override
210  {
211  KRATOS_TRY
212 
213  rComponentVariables.push_back(DISPLACEMENT);
214 
215  KRATOS_CATCH(" ")
216  }
217 
221 
225 
226  double& GetValue(const Variable<double> & , double& rValue) override;
227 
231 
233  virtual std::string Info() const override
234  {
235  std::stringstream buffer;
236  buffer << "SmallStrainUmatModel";
237  return buffer.str();
238  }
239 
241  virtual void PrintInfo(std::ostream& rOStream) const override
242  {
243  rOStream << "SmallStrainUmatModel";
244  }
245 
247  virtual void PrintData(std::ostream& rOStream) const override
248  {
249  rOStream << "SmallStrainUmatModel Data";
250  }
251 
252 
256 
257 
259 
260  protected:
263 
264 
268 
270 
274 
275 
279 
280 
284 
285 
286  virtual void InitializeStateVariables( Vector & rStateVariables, const Properties & rMaterialProperties)
287  {
288  };
289 
290 
291  //************//
292 
293  void InitializeElasticData(ModelDataType& rValues, UmatDataType& rVariables);
294 
295  /*
296  Get the dimension of StateVariables
297  */
298 
299  virtual unsigned int GetNumberOfStateVariables() {
300  KRATOS_ERROR << " Calling the base GetNumberOfStateVariables " << std::endl;
301  };
302 
303  /*
304  Create the vector with constitutive parameters value
305  */
306  virtual void CreateConstitutiveParametersVector(double* & rpVector, int & rNumberParameters, const Properties & rMaterialProperties) {
307  KRATOS_ERROR << " Calling the base CreateConstitutiveParametersVector " << std::endl;
308  };
309 
310  /*
311  Create state variables vector
312  ( tensor like should be "rotated" in large strains )
313  */
314  virtual void CreateStateVariablesVector(double * & rpStateVariables,int & rNumberStateVariables)
315  {
316  rpStateVariables = new double[rNumberStateVariables];
317  for (int i = 0; i < rNumberStateVariables; i++)
318  {
319  rpStateVariables[i] = mStateVariablesFinalized(i);
320  }
321  };
322 
323 
324 
325  /*
326  Create strain_n and incremental strain
327  ( for large strains should be overrided )
328  */
329 
330  virtual void CreateStrainsVectors( UmatDataType & rVariables, double* & rpStrain, double* & rpIncrementalStrain)
331  {
332 
333  VectorType StrainVector;
335 
336  rpStrain = new double[6];
337  rpIncrementalStrain = new double[6];
338 
339  for (unsigned int i = 0; i < 6; i++) {
340  rpStrain[i] = mStrainVectorFinalized(i);
341  rpIncrementalStrain[i] = ( StrainVector(i) - mStrainVectorFinalized(i) );
342 
343  }
344 
345  }
346 
347  /*
348  Create stress_n
349  ( for large strains should be overrided )
350  */
351  virtual void CreateStressAtInitialState( UmatDataType & rVariables, double* & rpStressVector)
352  {
353 
354  rpStressVector = new double[6];
355  for (unsigned int i = 0; i < 6; i++) {
356  rpStressVector[i] = mStressVectorFinalized(i);
357  }
358 
359  }
360 
361  /*
362  Update constitutive model variables
363  */
364  virtual void UpdateVariables( UmatDataType & rVariables, double* & rpStressVector, double* & rpStateVariables, double Pressure = 0.0)
365  {
366  VectorType StrainVector;
368  mStrainVectorFinalized = StrainVector;
369 
370 
371  if ( fabs( Pressure) < 1e-5) {
372  for (unsigned int i = 0; i < 6; i++)
373  mStressVectorFinalized(i) = rpStressVector[i];
374  } else {
375 
376  double meanP = 0;
377  for (unsigned int i = 0; i < 3; i++)
378  meanP += rpStressVector[i];
379  meanP /= 3.0;
380 
381  for (unsigned int i = 0; i < 3; i++)
382  mStressVectorFinalized(i) = rpStressVector[i] + Pressure - meanP;
383  for (unsigned int i = 3; i < 6; i++)
384  mStressVectorFinalized(i) = rpStressVector[i];
385  }
386 
387  int nStateVariables = this->GetNumberOfStateVariables();
388  for (int i = 0; i < nStateVariables; i++)
389  mStateVariablesFinalized(i) = rpStateVariables[i];
390 
391  }
392 
393  /*
394  Number of the constitutive equation in the fortran wrapper
395  */
397  {
398  KRATOS_ERROR << " Calling the base case of GetConstitutiveEquationNumber " << std::endl;
399  return 0;
400  }
401 
402  virtual void SetConstitutiveMatrix( Matrix & rC, const Matrix & rCBig, const MatrixType& rStressMatrix);
406 
407 
411 
412 
416 
417 
419 
420  private:
421 
424 
425 
429 
430 
434 
435 
439 
440 
444 
445 
449  friend class Serializer;
450 
451  virtual void save(Serializer& rSerializer) const override
452  {
454  rSerializer.save("InitializedModel", mInitializedModel );
455  rSerializer.save("StressVectorFinalized", mStressVectorFinalized );
456  rSerializer.save("StrainVectorFinalized", mStrainVectorFinalized );
457  rSerializer.save("StateVariablesFinalized", mStateVariablesFinalized );
458  }
459 
460  virtual void load(Serializer& rSerializer) override
461  {
463  rSerializer.load("InitializedModel", mInitializedModel );
464  rSerializer.load("StressVectorFinalized", mStressVectorFinalized );
465  rSerializer.load("StrainVectorFinalized", mStrainVectorFinalized );
466  rSerializer.load("StateVariablesFinalized", mStateVariablesFinalized );
467  }
468 
472 
473 
477 
479 
480  }; // Class SmallStrainUmatModel
481 
483 
486 
487 
491 
493 
495 
496 } // namespace Kratos.
497 
498 #endif // KRATOS_SMALL_STRAIN_UMAT_MODEL_H_INCLUDED defined
499 
500 
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
Short class definition.
Definition: constitutive_model.hpp:52
static void StrainTensorToVector(const MatrixType &rMatrix, array_1d< double, 6 > &rVector)
Definition: constitutive_model_utilities.hpp:646
Definition: flags.h:58
Definition: amatrix_interface.h:41
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_model.hpp:51
virtual void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: small_strain_umat_model.hpp:247
virtual void InitializeStateVariables(Vector &rStateVariables, const Properties &rMaterialProperties)
Definition: small_strain_umat_model.hpp:286
VectorType mStressVectorFinalized
Definition: small_strain_umat_model.hpp:272
virtual void SetValue(const Variable< Matrix > &rThisVariable, const Matrix &rValue, const ProcessInfo &rCurrentProcessInfo) override
Definition: small_strain_umat_model.hpp:187
virtual void CreateStateVariablesVector(double *&rpStateVariables, int &rNumberStateVariables)
Definition: small_strain_umat_model.hpp:314
virtual void CreateStressAtInitialState(UmatDataType &rVariables, double *&rpStressVector)
Definition: small_strain_umat_model.hpp:351
virtual std::string Info() const override
Turn back information as a string.
Definition: small_strain_umat_model.hpp:233
virtual void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: small_strain_umat_model.hpp:241
UmatModelData UmatDataType
Definition: small_strain_umat_model.hpp:87
virtual int GetConstitutiveEquationNumber()
Definition: small_strain_umat_model.hpp:396
virtual unsigned int GetNumberOfStateVariables()
Definition: small_strain_umat_model.hpp:299
virtual void GetDomainVariablesList(std::vector< Variable< double > > &rScalarVariables, std::vector< Variable< array_1d< double, 3 > > > &rComponentVariables) override
Definition: small_strain_umat_model.hpp:208
virtual void UpdateVariables(UmatDataType &rVariables, double *&rpStressVector, double *&rpStateVariables, double Pressure=0.0)
Definition: small_strain_umat_model.hpp:364
VectorType mStrainVectorFinalized
Definition: small_strain_umat_model.hpp:273
virtual void SetValue(const Variable< Vector > &rThisVariable, const Vector &rValue, const ProcessInfo &rCurrentProcessInfo) override
Definition: small_strain_umat_model.hpp:170
Vector mStateVariablesFinalized
Definition: small_strain_umat_model.hpp:271
bool mInitializedModel
Definition: small_strain_umat_model.hpp:269
virtual void CreateConstitutiveParametersVector(double *&rpVector, int &rNumberParameters, const Properties &rMaterialProperties)
Definition: small_strain_umat_model.hpp:306
virtual void CreateStrainsVectors(UmatDataType &rVariables, double *&rpStrain, double *&rpIncrementalStrain)
Definition: small_strain_umat_model.hpp:330
KRATOS_CLASS_POINTER_DEFINITION(SmallStrainUmatModel)
Pointer definition of SmallStrainUmatModel.
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_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
#define KRATOS_ERROR
Definition: exception.h:161
Matrix MatrixType
Definition: geometrical_transformation_utilities.h:55
double CalculateStrainEnergy(Element &rElement)
Definition: mpm_energy_calculation_utility.cpp:89
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
integer i
Definition: TensorModule.f:17
e
Definition: run_cpp_mpi_tests.py:31
Definition: constitutive_model_data.hpp:92
Definition: constitutive_model_data.hpp:383
const MaterialData & GetMaterialParameters() const
Definition: constitutive_model_data.hpp:462
Definition: small_strain_umat_model.hpp:55
const ModelDataType & GetModelData() const
Definition: small_strain_umat_model.hpp:71
const MaterialDataType & GetMaterialParameters() const
Definition: small_strain_umat_model.hpp:72
void SetModelData(const ModelDataType &rModelData)
Definition: small_strain_umat_model.hpp:68
const Flags & GetState() const
Definition: small_strain_umat_model.hpp:78
MatrixType TotalStrainMatrix
Definition: small_strain_umat_model.hpp:64
void SetState(Flags &rState)
Definition: small_strain_umat_model.hpp:67
MatrixType IncrementalDeformation
Definition: small_strain_umat_model.hpp:63
Flags & State()
Definition: small_strain_umat_model.hpp:75