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.
large_strain_3D_law.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosConstitutiveModelsApplication $
3 // Created by: $Author: JMCarbonell $
4 // Last modified by: $Co-Author: $
5 // Date: $Date: April 2017 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined (KRATOS_LARGE_STRAIN_3D_LAW_H_INCLUDED)
11 #define KRATOS_LARGE_STRAIN_3D_LAW_H_INCLUDED
12 
13 // System includes
14 
15 // External includes
16 
17 // Project includes
20 
21 namespace Kratos
22 {
27  class KRATOS_API(CONSTITUTIVE_MODELS_APPLICATION) LargeStrain3DLaw : public Constitutive3DLaw
28  {
29  public:
30 
33 
34  typedef ConstitutiveModel ModelType; //large_strain model
35  typedef ModelType::Pointer ModelTypePointer;
36 
39 
43 
46 
49 
51  LargeStrain3DLaw(const LargeStrain3DLaw& rOther);
52 
55 
57  ConstitutiveLaw::Pointer Clone() const override;
58 
60  ~LargeStrain3DLaw() override;
61 
65 
66 
70 
74  void InitializeMaterial(const Properties& rProperties,
75  const GeometryType& rElementGeometry,
76  const Vector& rShapeFunctionsValues ) override;
77 
84  void CalculateMaterialResponsePK2(Parameters & rValues) override;
85 
92  void CalculateMaterialResponseKirchhoff (Parameters & rValues) override;
93 
94 
99  void GetLawFeatures(Features& rFeatures) override;
100 
101 
106  void GetModelFeatures(Features& rFeatures);
107 
108 
118  int Check(const Properties& rProperties, const GeometryType& rElementGeometry, const ProcessInfo& rCurrentProcessInfo) const override;
119 
123 
128  bool Has( const Variable<double>& rThisVariable ) override;
129 
134  void SetValue( const Variable<double>& rVariable,
135  const double& rValue,
136  const ProcessInfo& rCurrentProcessInfo ) override;
137 
138  void SetValue( const Variable<Vector>& rVariable,
139  const Vector& rValue,
140  const ProcessInfo& rCurrentProcessInfo ) override;
141 
142  void SetValue( const Variable<Matrix>& rVariable,
143  const Matrix& rValue,
144  const ProcessInfo& rCurrentProcessInfo ) override;
145 
150  double& GetValue( const Variable<double>& rThisVariable, double& rValue ) override;
151 
152  Matrix& GetValue( const Variable<Matrix>& rThisVariable, Matrix& rValue ) override;
156 
160 
162  std::string Info() const override
163  {
164  std::stringstream buffer;
165  buffer << "LargeStrain3DLaw";
166  return buffer.str();
167  }
168 
170  void PrintInfo(std::ostream& rOStream) const override
171  {
172  rOStream << "LargeStrain3DLaw";
173  }
174 
176  void PrintData(std::ostream& rOStream) const override
177  {
178  rOStream << "LargeStrain3DLaw Data";
179  mpModel->PrintData(rOStream);
180  }
181 
182 
186 
187 
189 
190 
191  protected:
192 
195 
199 
200  //constitutive model
202 
203  //internal elastic variables
204 
205  //stored total deformation gradient for incremental strain update
208 
212 
216 
217 
225  void CalculateMaterialResponsePK2(Parameters & rValues, ModelDataType& rModelValues) override;
226 
234  void CalculateMaterialResponseKirchhoff (Parameters & rValues, ModelDataType& rModelValues) override;
235 
236 
240  void InitializeModelData(Parameters& rValues, ModelDataType& rModelValues) override;
241 
242 
246  void FinalizeModelData(Parameters& rValues, ModelDataType& rModelValues) override;
247 
253  virtual void CalculateStressVector(ModelDataType& rModelValues, Vector& rStressVector);
254 
255 
261  virtual void CalculateConstitutiveMatrix(ModelDataType& rModelValues, Matrix& rConstitutiveMatrix);
262 
263 
270  virtual void CalculateStressVectorAndConstitutiveMatrix(ModelDataType& rModelValues, Vector& rStressVector, Matrix& rConstitutiveMatrix);
271 
272 
274 
275  private:
276 
282 
286 
290 
295 
299  friend class Serializer;
300 
301  void save(Serializer& rSerializer) const override
302  {
304 
305  rSerializer.save("mpModel",mpModel);
306  rSerializer.save("mTotalDeformationDet",mTotalDeformationDet);
307  rSerializer.save("mInverseTotalDeformationMatrix",mInverseTotalDeformationMatrix);
308  }
309 
310  void load(Serializer& rSerializer) override
311  {
313 
314  rSerializer.load("mpModel",mpModel);
315  rSerializer.load("mTotalDeformationDet",mTotalDeformationDet);
316  rSerializer.load("mInverseTotalDeformationMatrix",mInverseTotalDeformationMatrix);
317  }
318 
319 
323 
324 
328 
329 
331  }; // Class LargeStrain3DLaw
332 
334 
337 
338 
342 
344 
346 
347 } // namespace Kratos.
348 #endif // KRATOS_LARGE_STRAIN_3D_LAW_H_INCLUDED defined
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
Definition: constitutive_3D_law.hpp:30
Short class definition.
Definition: constitutive_model.hpp:52
Geometry base class.
Definition: geometry.h:71
Definition: large_strain_3D_law.hpp:28
ModelTypePointer mpModel
Definition: large_strain_3D_law.hpp:201
ConstitutiveModel ModelType
Definition: large_strain_3D_law.hpp:34
KRATOS_CLASS_POINTER_DEFINITION(LargeStrain3DLaw)
Pointer definition of LargeStrain3DLaw.
double mTotalDeformationDet
Definition: large_strain_3D_law.hpp:206
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: large_strain_3D_law.hpp:170
std::string Info() const override
Turn back information as a string.
Definition: large_strain_3D_law.hpp:162
MatrixType mInverseTotalDeformationMatrix
Definition: large_strain_3D_law.hpp:207
ModelType::Pointer ModelTypePointer
Definition: large_strain_3D_law.hpp:35
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: large_strain_3D_law.hpp:176
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
#define KRATOS_SERIALIZE_SAVE_BASE_CLASS(Serializer, BaseType)
Definition: define.h:812
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
bool Has(const std::string &ModelerName)
Checks if the modeler is registered.
Definition: modeler_factory.cpp:24
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 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
Definition: constitutive_model_data.hpp:383