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.
rigid_body_element.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosContactMechanicsApplication $
3 // Created by: $Author: JMCarbonell $
4 // Last modified by: $Co-Author: $
5 // Date: $Date: July 2016 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_RIGID_BODY_ELEMENT_H_INCLUDED )
11 #define KRATOS_RIGID_BODY_ELEMENT_H_INCLUDED
12 
13 // System includes
14 
15 // External includes
16 
17 // Project includes
18 #include "includes/element.h"
20 
21 
22 namespace Kratos
23 {
38 
40 
46 class KRATOS_API(CONTACT_MECHANICS_APPLICATION) RigidBodyElement
47  :public Element
48 {
49 public:
50 
58  typedef Node NodeType;
65 
68 
70 
71 protected:
72 
76  KRATOS_DEFINE_LOCAL_FLAG( COMPUTE_RHS_VECTOR );
77  KRATOS_DEFINE_LOCAL_FLAG( COMPUTE_LHS_MATRIX );
78 
83  {
84  double Mass; // Mass of the Rigid Body
85  Matrix InertiaTensor; // Global Inertia Tensor
86  };
87 
92  {
93  private:
94 
95  //variables including all integration points
96  const ProcessInfo* pProcessInfo;
97 
98  public:
99 
100  //section properties
104 
105  void SetProcessInfo(const ProcessInfo& rProcessInfo)
106  {
107  pProcessInfo=&rProcessInfo;
108  }
109 
111  {
112  return *pProcessInfo;
113  }
114 
115  void Initialize(const unsigned int& dimension, const ProcessInfo& rProcessInfo)
116  {
117  noalias(VolumeForce) = ZeroVector(3);
118  DeltaPosition.resize(1,dimension,false);
119  noalias(DeltaPosition) = ZeroMatrix(1, dimension);
120  pProcessInfo=&rProcessInfo;
121  }
122 
123  };
124 
125 
134  {
135  private:
136 
137  //for calculation local system with compacted LHS and RHS
138  MatrixType *mpLeftHandSideMatrix;
139  VectorType *mpRightHandSideVector;
140 
141  public:
142 
143  //calculation flags
145 
149  void SetLeftHandSideMatrix( MatrixType& rLeftHandSideMatrix ) { mpLeftHandSideMatrix = &rLeftHandSideMatrix; };
150 
151  void SetRightHandSideVector( VectorType& rRightHandSideVector ) { mpRightHandSideVector = &rRightHandSideVector; };
152 
153 
157  MatrixType& GetLeftHandSideMatrix() { return *mpLeftHandSideMatrix; };
158 
159  VectorType& GetRightHandSideVector() { return *mpRightHandSideVector; };
160 
161  };
162 
163 
164 public:
165 
168 
171 
173  RigidBodyElement(IndexType NewId, GeometryType::Pointer pGeometry);
174 
175  RigidBodyElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties);
176 
177  RigidBodyElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties, NodesContainerType::Pointer pNodes);
178 
180  RigidBodyElement(RigidBodyElement const& rOther);
181 
183  virtual ~RigidBodyElement();
184 
185 
189 
190 
198  Element::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes, PropertiesType::Pointer pProperties) const override;
199 
207  Element::Pointer Clone(IndexType NewId, NodesArrayType const& ThisNodes) const override;
208 
209 
210 
211  //************* GETTING METHODS
212 
213 
217  void GetDofList(DofsVectorType& rElementalDofList, const ProcessInfo& rCurrentProcessInfo) const override;
218 
222  void EquationIdVector(EquationIdVectorType& rResult, const ProcessInfo& rCurrentProcessInfo) const override;
223 
227  void GetValuesVector(Vector& rValues, int Step = 0) const override;
228 
232  void GetFirstDerivativesVector(Vector& rValues, int Step = 0) const override;
233 
237  void GetSecondDerivativesVector(Vector& rValues, int Step = 0) const override;
238 
239 
240  //************* STARTING - ENDING METHODS
241 
246  void Initialize(const ProcessInfo& rCurrentProcessInfo) override;
247 
251  void InitializeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
252 
256  void InitializeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo) override;
257 
261  void FinalizeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo) override;
262 
263 
267  void FinalizeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
268 
269  //************* COMPUTING METHODS
270 
279  void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix, VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo) override;
280 
287  void CalculateRightHandSide(VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo) override;
288 
289 
296  void CalculateLeftHandSide(MatrixType& rLeftHandSideMatrix, const ProcessInfo& rCurrentProcessInfo) override;
297 
298 
306  void CalculateSecondDerivativesContributions(MatrixType& rLeftHandSideMatrix,
307  VectorType& rRightHandSideVector,
308  const ProcessInfo& rCurrentProcessInfo) override;
309 
316  void CalculateSecondDerivativesLHS(MatrixType& rLeftHandSideMatrix,
317  const ProcessInfo& rCurrentProcessInfo) override;
318 
319 
326  void CalculateSecondDerivativesRHS(VectorType& rRightHandSideVector,
327  const ProcessInfo& rCurrentProcessInfo) override;
328 
335  void CalculateMassMatrix(MatrixType& rMassMatrix, const ProcessInfo& rCurrentProcessInfo) override;
336 
337 
351  void AddExplicitContribution(const VectorType& rRHSVector, const Variable<VectorType>& rRHSVariable, const Variable<array_1d<double,3> >& rDestinationVariable, const ProcessInfo& rCurrentProcessInfo) override;
352 
360  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
361 
365 
369 
374  std::string Info() const override
375  {
376  std::stringstream buffer;
377  buffer << "Rigid Body Element #" << Id();
378  return buffer.str();
379  }
380 
382  void PrintInfo(std::ostream& rOStream) const override
383  {
384  rOStream << "Rigid Body Element #" << Id();
385  }
386 
388  void PrintData(std::ostream& rOStream) const override
389  {
390  GetGeometry().PrintData(rOStream);
391  }
396 
397 protected:
398 
404 
409 
410  NodesContainerType::Pointer mpNodes;
411 
415 
416 
420 
421 
425  void CalculateDynamicSystem(LocalSystemComponents& rLocalSystem,
426  const ProcessInfo& rCurrentProcessInfo);
427 
431  virtual void InitializeSystemMatrices(MatrixType& rLeftHandSideMatrix,
432  VectorType& rRightHandSideVector,
433  Flags& rCalculationFlags);
434 
438  ArrayType& MapToInitialLocalFrame(ArrayType& rVariable);
439 
440 
444  void CalculateRigidBodyProperties(RigidBodyProperties & rRigidBody);
445 
446 
450  virtual void CalculateAndAddLHS(MatrixType& rLeftHandSideMatrix,
451  ElementVariables& rVariables);
452 
456  virtual void CalculateAndAddRHS(VectorType& rRightHandSideVector,
457  ElementVariables& rVariables);
458 
462  virtual void CalculateAndAddExternalForces(VectorType& rRightHandSideVector,
463  ElementVariables& rVariables);
464 
468  virtual void CalculateAndAddInertiaLHS(MatrixType& rLeftHandSideMatrix,
469  ElementVariables& rVariables);
470 
474  virtual void CalculateAndAddInertiaRHS(VectorType& rRightHandSideVector,
475  ElementVariables& rVariables);
476 
477 
481  virtual void GetTimeIntegrationParameters(double& rP0,double& rP1,double& rP2,
482  const ProcessInfo& rCurrentProcessInfo);
483 
487  virtual ArrayType& CalculateVolumeForce(ArrayType& rVolumeForce);
488 
489 
493  virtual void CalculateRotationLinearPartTensor(ArrayType& rRotationVector, Matrix& rRotationTensor);
494 
495 
499  virtual void UpdateRigidBodyNodes(const ProcessInfo& rCurrentProcessInfo);
500 
504  virtual SizeType GetDofsSize() const;
505 
509  virtual void MapLocalToGlobalSystem(LocalSystemComponents& rLocalSystem);
510 
521 
522 private:
523 
542  friend class Serializer;
543 
544  virtual void save(Serializer& rSerializer) const override;
545 
546  virtual void load(Serializer& rSerializer) override;
547 
554 
555 
556 }; // Class RigidBodyElement
557 
558 } // namespace Kratos.
559 #endif // KRATOS_RIGID_BODY_ELEMENT_H_INCLUDED defined
Definition: beam_math_utilities.hpp:31
Base class for all Elements.
Definition: element.h:60
std::size_t SizeType
Definition: element.h:94
Definition: flags.h:58
std::size_t SizeType
Definition: geometry_data.h:173
void resize(std::size_t NewSize1, std::size_t NewSize2, bool preserve=0)
Definition: amatrix_interface.h:224
This class defines the node.
Definition: node.h:65
A sorted associative container similar to an STL set, but uses a vector to store pointers to its data...
Definition: pointer_vector_set.h:72
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
Rigid Body Element for 3D space dimension.
Definition: rigid_body_element.hpp:48
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: rigid_body_element.hpp:382
KRATOS_DEFINE_LOCAL_FLAG(COMPUTE_LHS_MATRIX)
std::string Info() const override
Turn back information as a string.
Definition: rigid_body_element.hpp:374
Quaternion< double > QuaternionType
Type definition for quaternion.
Definition: rigid_body_element.hpp:56
GeometryData::SizeType SizeType
Type for size.
Definition: rigid_body_element.hpp:62
BeamMathUtils< double > BeamMathUtilsType
Definition: rigid_body_element.hpp:54
QuaternionType mInitialLocalQuaternion
Definition: rigid_body_element.hpp:408
array_1d< double, 3 > ArrayType
Type of vector.
Definition: rigid_body_element.hpp:64
KRATOS_DEFINE_LOCAL_FLAG(COMPUTE_RHS_VECTOR)
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: rigid_body_element.hpp:388
PointerVectorSet< NodeType, IndexedObject > NodesContainerType
Type for nodes container.
Definition: rigid_body_element.hpp:60
Node NodeType
Type for nodes.
Definition: rigid_body_element.hpp:58
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(RigidBodyElement)
Counted pointer of RigidBodyElement.
NodesContainerType::Pointer mpNodes
Definition: rigid_body_element.hpp:410
RigidBodyElement()
Serializer constructor.
Definition: rigid_body_element.hpp:170
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
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 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
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
KratosZeroVector< double > ZeroVector
Definition: amatrix_interface.h:561
KratosZeroMatrix< double > ZeroMatrix
Definition: amatrix_interface.h:559
ModelPart::NodesContainerType NodesArrayType
Definition: gid_gauss_point_container.h:42
T & noalias(T &TheMatrix)
Definition: amatrix_interface.h:484
int dimension
Definition: isotropic_damage_automatic_differentiation.py:123
def load(f)
Definition: ode_solve.py:307
Definition: rigid_body_element.hpp:92
const ProcessInfo & GetProcessInfo()
Definition: rigid_body_element.hpp:110
Matrix DeltaPosition
Definition: rigid_body_element.hpp:103
ArrayType VolumeForce
Definition: rigid_body_element.hpp:102
RigidBodyProperties RigidBody
Definition: rigid_body_element.hpp:101
void SetProcessInfo(const ProcessInfo &rProcessInfo)
Definition: rigid_body_element.hpp:105
void Initialize(const unsigned int &dimension, const ProcessInfo &rProcessInfo)
Definition: rigid_body_element.hpp:115
Definition: rigid_body_element.hpp:134
Flags CalculationFlags
Definition: rigid_body_element.hpp:144
void SetRightHandSideVector(VectorType &rRightHandSideVector)
Definition: rigid_body_element.hpp:151
void SetLeftHandSideMatrix(MatrixType &rLeftHandSideMatrix)
Definition: rigid_body_element.hpp:149
MatrixType & GetLeftHandSideMatrix()
Definition: rigid_body_element.hpp:157
VectorType & GetRightHandSideVector()
Definition: rigid_body_element.hpp:159
Definition: rigid_body_element.hpp:83
Matrix InertiaTensor
Definition: rigid_body_element.hpp:85
double Mass
Definition: rigid_body_element.hpp:84