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.
updated_lagrangian.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: Vicente Mataix Ferrandiz
10 //
11 
12 #pragma once
13 
14 // System includes
15 
16 // External includes
17 
18 // Project includes
19 #include "includes/define.h"
20 #include "custom_elements/base_solid_element.h"
21 #include "includes/variables.h"
23 
24 namespace Kratos
25 {
34 
38 
42 
51 class KRATOS_API(STRUCTURAL_MECHANICS_APPLICATION) UpdatedLagrangian
52  : public BaseSolidElement
53 {
54 public:
57 
60 
62  typedef ConstitutiveLawType::Pointer ConstitutiveLawPointerType;
63 
66 
69 
71  typedef std::size_t IndexType;
72 
74  typedef std::size_t SizeType;
75 
78 
82 
84  UpdatedLagrangian(IndexType NewId, GeometryType::Pointer pGeometry);
85  UpdatedLagrangian(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties);
86 
87  // Copy constructor
89  :BaseType(rOther)
90  ,mF0Computed(rOther.mF0Computed)
91  ,mDetF0(rOther.mDetF0)
92  ,mF0(rOther.mF0)
93  {};
94 
96  ~UpdatedLagrangian() override;
97 
104 
109  void Initialize(const ProcessInfo& rCurrentProcessInfo) override;
110 
115  void InitializeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
116 
121  void FinalizeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
122 
130  Element::Pointer Clone (
131  IndexType NewId,
132  NodesArrayType const& rThisNodes
133  ) const override;
134 
142  Element::Pointer Create(
143  IndexType NewId,
144  GeometryType::Pointer pGeom,
145  PropertiesType::Pointer pProperties
146  ) const override;
147 
155  Element::Pointer Create(
156  IndexType NewId,
157  NodesArrayType const& ThisNodes,
158  PropertiesType::Pointer pProperties
159  ) const override;
160 
161 
169  const Variable<bool>& rVariable,
170  std::vector<bool>& rOutput,
171  const ProcessInfo& rCurrentProcessInfo
172  ) override;
173 
181  const Variable<int>& rVariable,
182  std::vector<int>& rOutput,
183  const ProcessInfo& rCurrentProcessInfo
184  ) override;
185 
193  const Variable<double>& rVariable,
194  std::vector<double>& rOutput,
195  const ProcessInfo& rCurrentProcessInfo
196  ) override;
197 
205  const Variable<array_1d<double, 3>>& rVariable,
206  std::vector<array_1d<double, 3>>& rOutput,
207  const ProcessInfo& rCurrentProcessInfo
208  ) override;
209 
217  const Variable<array_1d<double, 6>>& rVariable,
218  std::vector<array_1d<double, 6>>& rOutput,
219  const ProcessInfo& rCurrentProcessInfo
220  ) override;
221 
229  const Variable<Vector>& rVariable,
230  std::vector<Vector>& rOutput,
231  const ProcessInfo& rCurrentProcessInfo
232  ) override;
233 
241  const Variable<Matrix>& rVariable,
242  std::vector<Matrix>& rOutput,
243  const ProcessInfo& rCurrentProcessInfo
244  ) override;
245 
253  const Variable<double>& rVariable,
254  const std::vector<double>& rValues,
255  const ProcessInfo& rCurrentProcessInfo
256  ) override;
257 
265  const Variable<Matrix>& rVariable,
266  const std::vector<Matrix>& rValues,
267  const ProcessInfo& rCurrentProcessInfo
268  ) override;
269 
273 
280 
282  std::string Info() const override
283  {
284  std::stringstream buffer;
285  buffer << "Updated Lagrangian Solid Element #" << Id() << "\nConstitutive law: " << BaseType::mConstitutiveLawVector[0]->Info();
286  return buffer.str();
287  }
288 
290  void PrintInfo(std::ostream& rOStream) const override
291  {
292  rOStream << "Updated Lagrangian Solid Element #" << Id() << "\nConstitutive law: " << BaseType::mConstitutiveLawVector[0]->Info();
293  }
294 
296  void PrintData(std::ostream& rOStream) const override
297  {
298  pGetGeometry()->PrintData(rOStream);
299  }
300 
305 
306 protected:
309 
313 
314  /* Historical total elastic deformation measure */
315  bool mF0Computed; // To avoid computing more than once the historical total elastic deformation measure
316  std::vector<double> mDetF0; // The historical total elastic deformation measure determinant
317  std::vector<Matrix> mF0; // The historical total elastic deformation measure
318 
322 
324  {
325  }
326 
334  const bool rF0Computed,
335  const std::vector<double>& rDetF0,
336  const std::vector<Matrix>& rF0
337  )
338  {
339  mF0Computed = rF0Computed;
340  mDetF0 = rDetF0;
341  mF0 = rF0;
342  }
343 
347  ConstitutiveLaw::StressMeasure GetStressMeasure() const override;
348 
354  void UpdateHistoricalDatabase(
355  KinematicVariables& rThisKinematicVariables,
356  const IndexType PointNumber
357  );
358 
367  void CalculateAll(
368  MatrixType& rLeftHandSideMatrix,
369  VectorType& rRightHandSideVector,
370  const ProcessInfo& rCurrentProcessInfo,
371  const bool CalculateStiffnessMatrixFlag,
372  const bool CalculateResidualVectorFlag
373  ) override;
374 
381  void CalculateKinematicVariables(
382  KinematicVariables& rThisKinematicVariables,
383  const IndexType PointNumber,
384  const GeometryType::IntegrationMethod& rIntegrationMethod
385  ) override;
386 
396  double CalculateDerivativesOnReferenceConfiguration(
397  Matrix& J0,
398  Matrix& InvJ0,
399  Matrix& DN_DX,
400  const IndexType PointNumber,
401  IntegrationMethod ThisIntegrationMethod
402  ) const override;
403 
417 
418 private:
421 
425 
429 
437  void CalculateB(
438  Matrix& rB,
439  const Matrix& rDN_DX,
440  const SizeType StrainSize,
441  const IndexType PointNumber
442  );
443 
449  double ReferenceConfigurationDeformationGradientDeterminant(const IndexType PointNumber) const;
450 
456  Matrix ReferenceConfigurationDeformationGradient(const IndexType PointNumber) const;
457 
461 
462 
467 
471  friend class Serializer;
472 
473  // A private default constructor necessary for serialization
474 
475  void save(Serializer& rSerializer) const override;
476 
477  void load(Serializer& rSerializer) override;
478 
485  //UpdatedLagrangian& operator=(const UpdatedLagrangian& rOther);
487  //UpdatedLagrangian(const UpdatedLagrangian& rOther);
489 
490 }; // Class UpdatedLagrangian
491 
499 
500 } // namespace Kratos.
This is base class used to define the solid elements.
Definition: base_solid_element.h:67
Definition: constitutive_law.h:47
StressMeasure
Definition: constitutive_law.h:69
std::size_t IndexType
Definition: flags.h:74
This defines the geometrical object, base definition of the element and condition entities.
Definition: geometrical_object.h:58
IntegrationMethod
Definition: geometry_data.h:76
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
Updated Lagrangian element for 2D and 3D geometries.
Definition: updated_lagrangian.h:53
UpdatedLagrangian()
Definition: updated_lagrangian.h:323
std::string Info() const override
Turn back information as a string.
Definition: updated_lagrangian.h:282
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: updated_lagrangian.h:290
UpdatedLagrangian(UpdatedLagrangian const &rOther)
Definition: updated_lagrangian.h:88
ConstitutiveLaw ConstitutiveLawType
Reference type definition for constitutive laws.
Definition: updated_lagrangian.h:59
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(UpdatedLagrangian)
Counted pointer of UpdatedLagrangian.
std::size_t IndexType
The definition of the index type.
Definition: updated_lagrangian.h:71
bool mF0Computed
Definition: updated_lagrangian.h:315
std::vector< Matrix > mF0
Definition: updated_lagrangian.h:317
BaseSolidElement BaseType
The base element type.
Definition: updated_lagrangian.h:68
std::vector< double > mDetF0
Definition: updated_lagrangian.h:316
GeometryData::IntegrationMethod IntegrationMethod
Type definition for integration methods.
Definition: updated_lagrangian.h:65
ConstitutiveLawType::Pointer ConstitutiveLawPointerType
Pointer type for constitutive laws.
Definition: updated_lagrangian.h:62
void CloneUpdatedLagrangianDatabase(const bool rF0Computed, const std::vector< double > &rDetF0, const std::vector< Matrix > &rF0)
This method clones the element database.
Definition: updated_lagrangian.h:333
std::size_t SizeType
The definition of the sizetype.
Definition: updated_lagrangian.h:74
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: updated_lagrangian.h:296
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
std::size_t IndexType
The definition of the index type.
Definition: key_hash.h:35
Matrix MatrixType
Definition: geometrical_transformation_utilities.h:55
Modeler::Pointer Create(const std::string &ModelerName, Model &rModel, const Parameters ModelParameters)
Checks if the modeler is registered.
Definition: modeler_factory.cpp:30
void SetValuesOnIntegrationPoints(TObject &dummy, const Variable< TDataType > &rVariable, const std::vector< TDataType > &values, const ProcessInfo &rCurrentProcessInfo)
Definition: add_mesh_to_python.cpp:185
pybind11::list CalculateOnIntegrationPoints(TObject &dummy, const Variable< TDataType > &rVariable, const ProcessInfo &rProcessInfo)
Definition: add_mesh_to_python.cpp:142
void InitializeSolutionStep(ConstructionUtility &rThisUtil, std::string ThermalSubModelPartName, std::string MechanicalSubModelPartName, std::string HeatFluxSubModelPartName, std::string HydraulicPressureSubModelPartName, bool thermal_conditions, bool mechanical_conditions, int phase)
Definition: add_custom_utilities_to_python.cpp:45
void CalculateB(const GeometricalObject &rElement, const TMatrixType1 &rDN_DX, TMatrixType2 &rB)
This method computes the deformation tensor B (for small deformation solid elements)
Definition: structural_mechanics_element_utilities.h:105
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
std::size_t SizeType
The definition of the size type.
Definition: mortar_classes.h:43
ModelPart::NodesContainerType NodesArrayType
Definition: gid_gauss_point_container.h:42
def load(f)
Definition: ode_solve.py:307