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.
mpm_updated_lagrangian_UP.hpp
Go to the documentation of this file.
1 // | / |
2 // ' / __| _` | __| _ \ __|
3 // . \ | ( | | ( |\__ `
4 // _|\_\_| \__,_|\__|\___/ ____/
5 // Multi-Physics
6 //
7 // License: BSD License
8 // Kratos default license: kratos/license.txt
9 //
10 // Main authors: Ilaria Iaconeta
11 //
12 
13 
14 #if !defined(KRATOS_UPDATED_LAGRANGIAN_UP_H_INCLUDED )
15 #define KRATOS_UPDATED_LAGRANGIAN_UP_H_INCLUDED
16 
17 // System includes
18 
19 // External includes
20 
21 // Project includes
23 
24 namespace Kratos
25 {
40 
42 
49  : public MPMUpdatedLagrangian
50 {
51 public:
52 
58  typedef ConstitutiveLawType::Pointer ConstitutiveLawPointerType;
63 
67 
68 
69 
70 
71 public:
72 
73 
76 
79 
80 
82  MPMUpdatedLagrangianUP(IndexType NewId, GeometryType::Pointer pGeometry);
83 
84  MPMUpdatedLagrangianUP(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties);
85 
88 
90  ~MPMUpdatedLagrangianUP() override;
91 
95 
98 
102 
110  Element::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes, PropertiesType::Pointer pProperties) const override;
111 
112  Element::Pointer Create(IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties) const override;
113 
121  Element::Pointer Clone(IndexType NewId, NodesArrayType const& ThisNodes) const override;
122 
123 
124  //************* GETTING METHODS
125 
129  void GetDofList(DofsVectorType& rElementalDofList, const ProcessInfo& rCurrentProcessInfo) const override;
130 
134  void EquationIdVector(EquationIdVectorType& rResult, const ProcessInfo& rCurrentProcessInfo) const override;
135 
139  void GetValuesVector(Vector& rValues, int Step = 0) const override;
140 
144  void GetFirstDerivativesVector(Vector& rValues, int Step = 0) const override;
145 
149  void GetSecondDerivativesVector(Vector& rValues, int Step = 0) const override;
150 
151  //************* STARTING - ENDING METHODS
152 
157  void Initialize(const ProcessInfo& rCurrentProcessInfo) override;
158 
162  void InitializeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
163 
164 
165  //************* COMPUTING METHODS
166 
173  void CalculateMassMatrix(MatrixType& rMassMatrix,
174  const ProcessInfo& rCurrentProcessInfo) override;
175 
176 
177  //************************************************************************************
178  //************************************************************************************
186  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
187 
188 
192 
193  void CalculateOnIntegrationPoints(const Variable<double>& rVariable,
194  std::vector<double>& rValues,
195  const ProcessInfo& rCurrentProcessInfo) override;
196 
198  const Variable<double>& rVariable,
199  const std::vector<double>& rValues,
200  const ProcessInfo& rCurrentProcessInfo) override;
201 
209  std::string Info() const override
210  {
211  std::stringstream buffer;
212  buffer << "MPM Element #" << Id();
213  return buffer.str();
214  }
215 
217  void PrintInfo(std::ostream& rOStream) const override
218  {
219  rOStream << "MPM Element #" << Id();
220  }
221 
223  void PrintData(std::ostream& rOStream) const override
224  {
225  GetGeometry().PrintData(rOStream);
226  }
231 
232 protected:
241 
243  return GetGeometry().WorkingSpaceDimension() + 1;
244  }
245 
251  //virtual void CalculateElementalSystem(LocalSystemComponents& rLocalSystem,
252  //ProcessInfo& rCurrentProcessInfo);
256  void FinalizeStepVariables(GeneralVariables & rVariables, const ProcessInfo& rCurrentProcessInfo) override;
257 
262  void CalculateAndAddLHS(MatrixType& rLeftHandSideMatrix,
263  GeneralVariables& rVariables,
264  const double& rIntegrationWeight,
265  const ProcessInfo& rCurrentProcessInfo) override;
266 
271  void CalculateAndAddRHS(VectorType& rRightHandSideVector,
272  GeneralVariables& rVariables,
273  Vector& rVolumeForce,
274  const double& rIntegrationWeight,
275  const ProcessInfo& rCurrentProcessInfo) override;
276 
277 
282  void CalculateAndAddKuum(MatrixType& rLeftHandSideMatrix,
283  GeneralVariables& rVariables,
284  const double& rIntegrationWeight) override;
285 
289  void CalculateAndAddKuugUP(MatrixType& rLeftHandSideMatrix,
290  GeneralVariables& rVariables,
291  const double& rIntegrationWeight);
292 
296  virtual void CalculateAndAddKup (MatrixType& rK,
297  GeneralVariables & rVariables,
298  const double& rIntegrationWeight
299  );
300 
304  virtual void CalculateAndAddKpu(MatrixType& rK,
305  GeneralVariables & rVariables,
306  const double& rIntegrationWeight
307  );
308 
312  virtual void CalculateAndAddKpp(MatrixType& rK,
313  GeneralVariables & rVariables,
314  const double& rIntegrationWeight
315  );
316 
320  virtual void CalculateAndAddKppStab(MatrixType& rK,
321  GeneralVariables & rVariables,
322  const double& rIntegrationWeight
323  );
327  void CalculateAndAddExternalForces(VectorType& rRightHandSideVector,
328  GeneralVariables& rVariables,
329  Vector& rVolumeForce,
330  const double& rIntegrationWeight) override;
331 
335  void CalculateAndAddInternalForces(VectorType& rRightHandSideVector,
336  GeneralVariables & rVariables,
337  const double& rIntegrationWeight) override;
338 
342  virtual void CalculateAndAddPressureForces(VectorType& rRightHandSideVector,
343  GeneralVariables & rVariables,
344  const double& rIntegrationWeight
345  );
346 
350  virtual void CalculateAndAddStabilizedPressure(VectorType& rRightHandSideVector,
351  GeneralVariables & rVariables,
352  const double& rIntegrationWeight
353  );
354 
358  void CalculateKinematics(GeneralVariables& rVariables, const ProcessInfo& rCurrentProcessInfo) override;
359 
363  void InitializeGeneralVariables(GeneralVariables & rVariables, const ProcessInfo& rCurrentProcessInfo) override;
368  void UpdateGaussPoint(GeneralVariables & rVariables, const ProcessInfo& rCurrentProcessInfo) override;
369 
370 
374  virtual double& CalculatePUCoefficient(double& rCoefficient, GeneralVariables & rVariables);
375 
379  virtual double& CalculatePUDeltaCoefficient(double& rCoefficient, GeneralVariables & rVariables);
380 
384  void GetHistoricalVariables( GeneralVariables& rVariables) override;
385 
391  Matrix& rF,
392  Matrix& rDN_DX);
393 
397  double& CalculateVolumeChange(double& rVolumeChange, GeneralVariables& rVariables) override;
398 
409 
410 private:
411 
417 
418  double m_mp_pressure;
419 
423 
424 
428 
429 
434 
438  friend class Serializer;
439 
440  // A private default constructor necessary for serialization
441 
442  void save(Serializer& rSerializer) const override;
443 
444  void load(Serializer& rSerializer) override;
445 
446 
453 
454 }; // Class MPMUpdatedLagrangian
455 
463 
464 } // namespace Kratos.
465 #endif // KRATOS_UPDATED_LAGRANGIAN_H_INCLUDED defined
Definition: constitutive_law.h:47
StressMeasure
Definition: constitutive_law.h:69
std::size_t SizeType
Definition: element.h:94
std::vector< DofType::Pointer > DofsVectorType
Definition: element.h:100
Matrix MatrixType
Definition: element.h:90
std::vector< std::size_t > EquationIdVectorType
Definition: element.h:98
std::size_t IndexType
Definition: flags.h:74
GeometryType & GetGeometry()
Returns the reference of the geometry.
Definition: geometrical_object.h:158
IntegrationMethod
Definition: geometry_data.h:76
virtual void PrintData(std::ostream &rOStream) const
Print object's data.
Definition: geometry.h:3834
SizeType WorkingSpaceDimension() const
Definition: geometry.h:1287
IndexType Id() const
Definition: indexed_object.h:107
Large Displacement Lagrangian Element for 3D and 2D geometries. (base class)
Definition: mpm_updated_lagrangian.hpp:55
virtual void CalculateDeformationMatrix(Matrix &rB, const Matrix &rDN_DX, const Matrix &rN, const bool IsAxisymmetric=false)
Definition: mpm_updated_lagrangian.cpp:399
Large Displacement Lagrangian Element for 3D and 2D geometries. (base class)
Definition: mpm_updated_lagrangian_UP.hpp:50
ConstitutiveLaw ConstitutiveLawType
Definition: mpm_updated_lagrangian_UP.hpp:56
virtual void CalculateAndAddKup(MatrixType &rK, GeneralVariables &rVariables, const double &rIntegrationWeight)
Definition: mpm_updated_lagrangian_UP.cpp:759
void CalculateAndAddLHS(MatrixType &rLeftHandSideMatrix, GeneralVariables &rVariables, const double &rIntegrationWeight, const ProcessInfo &rCurrentProcessInfo) override
Definition: mpm_updated_lagrangian_UP.cpp:639
int Check(const ProcessInfo &rCurrentProcessInfo) const override
Definition: mpm_updated_lagrangian_UP.cpp:1220
void CalculateAndAddRHS(VectorType &rRightHandSideVector, GeneralVariables &rVariables, Vector &rVolumeForce, const double &rIntegrationWeight, const ProcessInfo &rCurrentProcessInfo) override
Definition: mpm_updated_lagrangian_UP.cpp:414
void GetHistoricalVariables(GeneralVariables &rVariables) override
Definition: mpm_updated_lagrangian_UP.cpp:1135
virtual void CalculateAndAddKppStab(MatrixType &rK, GeneralVariables &rVariables, const double &rIntegrationWeight)
Definition: mpm_updated_lagrangian_UP.cpp:865
MPMUpdatedLagrangianUP()
Empty constructor needed for serialization.
Definition: mpm_updated_lagrangian_UP.cpp:42
void CalculateKinematics(GeneralVariables &rVariables, const ProcessInfo &rCurrentProcessInfo) override
Definition: mpm_updated_lagrangian_UP.cpp:242
void GetFirstDerivativesVector(Vector &rValues, int Step=0) const override
Definition: mpm_updated_lagrangian_UP.cpp:1075
void CalculateAndAddInternalForces(VectorType &rRightHandSideVector, GeneralVariables &rVariables, const double &rIntegrationWeight) override
Definition: mpm_updated_lagrangian_UP.cpp:471
MPMUpdatedLagrangianUP & operator=(MPMUpdatedLagrangianUP const &rOther)
Assignment operator.
Definition: mpm_updated_lagrangian_UP.cpp:82
double & CalculateVolumeChange(double &rVolumeChange, GeneralVariables &rVariables) override
Definition: mpm_updated_lagrangian_UP.cpp:929
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: mpm_updated_lagrangian_UP.hpp:217
GeometryData::IntegrationMethod IntegrationMethod
Type definition for integration methods.
Definition: mpm_updated_lagrangian_UP.hpp:62
virtual void CalculateAndAddStabilizedPressure(VectorType &rRightHandSideVector, GeneralVariables &rVariables, const double &rIntegrationWeight)
Definition: mpm_updated_lagrangian_UP.cpp:576
void InitializeSolutionStep(const ProcessInfo &rCurrentProcessInfo) override
Definition: mpm_updated_lagrangian_UP.cpp:346
SizeType GetNumberOfDofs() override
Definition: mpm_updated_lagrangian_UP.hpp:242
Element::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
Definition: mpm_updated_lagrangian_UP.cpp:93
ConstitutiveLawType::StressMeasure StressMeasureType
StressMeasure from constitutive laws.
Definition: mpm_updated_lagrangian_UP.hpp:60
std::string Info() const override
Turn back information as a string.
Definition: mpm_updated_lagrangian_UP.hpp:209
void CalculateAndAddKuugUP(MatrixType &rLeftHandSideMatrix, GeneralVariables &rVariables, const double &rIntegrationWeight)
Definition: mpm_updated_lagrangian_UP.cpp:717
ConstitutiveLawType::Pointer ConstitutiveLawPointerType
Pointer type for constitutive laws.
Definition: mpm_updated_lagrangian_UP.hpp:58
virtual double & CalculatePUCoefficient(double &rCoefficient, GeneralVariables &rVariables)
Definition: mpm_updated_lagrangian_UP.cpp:498
virtual void CalculateAndAddPressureForces(VectorType &rRightHandSideVector, GeneralVariables &rVariables, const double &rIntegrationWeight)
Definition: mpm_updated_lagrangian_UP.cpp:531
virtual void CalculateAndAddKpu(MatrixType &rK, GeneralVariables &rVariables, const double &rIntegrationWeight)
Definition: mpm_updated_lagrangian_UP.cpp:790
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: mpm_updated_lagrangian_UP.hpp:223
void CalculateAndAddKuum(MatrixType &rLeftHandSideMatrix, GeneralVariables &rVariables, const double &rIntegrationWeight) override
Definition: mpm_updated_lagrangian_UP.cpp:678
void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const override
Definition: mpm_updated_lagrangian_UP.cpp:943
void SetValuesOnIntegrationPoints(const Variable< double > &rVariable, const std::vector< double > &rValues, const ProcessInfo &rCurrentProcessInfo) override
Definition: mpm_updated_lagrangian_UP.cpp:1193
void CalculateAndAddExternalForces(VectorType &rRightHandSideVector, GeneralVariables &rVariables, Vector &rVolumeForce, const double &rIntegrationWeight) override
Definition: mpm_updated_lagrangian_UP.cpp:445
~MPMUpdatedLagrangianUP() override
Destructor.
Definition: mpm_updated_lagrangian_UP.cpp:123
void CalculateOnIntegrationPoints(const Variable< double > &rVariable, std::vector< double > &rValues, const ProcessInfo &rCurrentProcessInfo) override
Definition: mpm_updated_lagrangian_UP.cpp:1177
virtual void CalculateDeformationMatrix(Matrix &rB, const Matrix &rDN_DX, const Matrix &rN, const bool IsAxisymmetric=false)
Definition: mpm_updated_lagrangian.cpp:399
Element::Pointer Clone(IndexType NewId, NodesArrayType const &ThisNodes) const override
Definition: mpm_updated_lagrangian_UP.cpp:106
void GetValuesVector(Vector &rValues, int Step=0) const override
Definition: mpm_updated_lagrangian_UP.cpp:1043
void Initialize(const ProcessInfo &rCurrentProcessInfo) override
Definition: mpm_updated_lagrangian_UP.cpp:131
void FinalizeStepVariables(GeneralVariables &rVariables, const ProcessInfo &rCurrentProcessInfo) override
Definition: mpm_updated_lagrangian_UP.cpp:1144
void UpdateGaussPoint(GeneralVariables &rVariables, const ProcessInfo &rCurrentProcessInfo) override
Definition: mpm_updated_lagrangian_UP.cpp:168
void CalculateMassMatrix(MatrixType &rMassMatrix, const ProcessInfo &rCurrentProcessInfo) override
Definition: mpm_updated_lagrangian_UP.cpp:997
void InitializeGeneralVariables(GeneralVariables &rVariables, const ProcessInfo &rCurrentProcessInfo) override
Definition: mpm_updated_lagrangian_UP.cpp:152
void GetSecondDerivativesVector(Vector &rValues, int Step=0) const override
Definition: mpm_updated_lagrangian_UP.cpp:1104
virtual double & CalculatePUDeltaCoefficient(double &rCoefficient, GeneralVariables &rVariables)
Definition: mpm_updated_lagrangian_UP.cpp:514
void GetDofList(DofsVectorType &rElementalDofList, const ProcessInfo &rCurrentProcessInfo) const override
Definition: mpm_updated_lagrangian_UP.cpp:974
virtual void CalculateAndAddKpp(MatrixType &rK, GeneralVariables &rVariables, const double &rIntegrationWeight)
Definition: mpm_updated_lagrangian_UP.cpp:821
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(MPMUpdatedLagrangianUP)
Counted pointer of LargeDisplacementElement.
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
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
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
def load(f)
Definition: ode_solve.py:307