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.
translatory_rigid_body_segregated_V_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: October 2018 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_TRANSLATORY_RIGID_BODY_SEGREGATED_V_ELEMENT_H_INCLUDED )
11 #define KRATOS_TRANSLATORY_RIGID_BODY_SEGREGATED_V_ELEMENT_H_INCLUDED
12 
13 // System includes
14 
15 // External includes
16 
17 // Project includes
19 
20 namespace Kratos
21 {
36 
38 
44 class KRATOS_API(CONTACT_MECHANICS_APPLICATION) TranslatoryRigidBodySegregatedVElement
46 {
47 public:
48 
56  typedef Node NodeType;
61 
64 
65  enum StepType{VELOCITY_STEP = 0, PRESSURE_STEP = 1};
66 
70 
73 
75  TranslatoryRigidBodySegregatedVElement(IndexType NewId, GeometryType::Pointer pGeometry);
76 
77  TranslatoryRigidBodySegregatedVElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties);
78 
79  TranslatoryRigidBodySegregatedVElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties, NodesContainerType::Pointer pNodes);
80 
83 
86 
90 
98  Element::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes, PropertiesType::Pointer pProperties) const override;
99 
107  Element::Pointer Clone(IndexType NewId, NodesArrayType const& ThisNodes) const override;
108 
109 
110 
111  //************* GETTING METHODS
112 
113 
117  void GetDofList(DofsVectorType& rElementalDofList, const ProcessInfo& rCurrentProcessInfo) const override;
118 
122  void EquationIdVector(EquationIdVectorType& rResult, const ProcessInfo& rCurrentProcessInfo) const override;
123 
127  void GetValuesVector(Vector& rValues, int Step = 0) const override;
128 
132  void GetFirstDerivativesVector(Vector& rValues, int Step = 0) const override;
133 
137  void GetSecondDerivativesVector(Vector& rValues, int Step = 0) const override;
138 
139 
140  //************* STARTING - ENDING METHODS
141 
146  void Initialize(const ProcessInfo& rCurrentProcessInfo) override;
147 
151  void InitializeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
152 
156  void InitializeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo) override;
157 
161  void FinalizeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo) override;
162 
163 
167  void FinalizeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
168 
169 
170  //************* COMPUTING METHODS
171 
180  void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix, VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo) override;
181 
188  void CalculateRightHandSide(VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo) override;
189 
190 
197  void CalculateLeftHandSide(MatrixType& rLeftHandSideMatrix, const ProcessInfo& rCurrentProcessInfo) override;
198 
199 
207  void CalculateSecondDerivativesContributions(MatrixType& rLeftHandSideMatrix,
208  VectorType& rRightHandSideVector,
209  const ProcessInfo& rCurrentProcessInfo) override;
210 
217  void CalculateSecondDerivativesLHS(MatrixType& rLeftHandSideMatrix,
218  const ProcessInfo& rCurrentProcessInfo) override;
219 
220 
227  void CalculateSecondDerivativesRHS(VectorType& rRightHandSideVector,
228  const ProcessInfo& rCurrentProcessInfo) override;
229 
236  void CalculateMassMatrix(MatrixType& rMassMatrix, const ProcessInfo& rCurrentProcessInfo) override;
237 
238  //************************************************************************************
239  //************************************************************************************
247  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
248 
252 
256 
261  std::string Info() const override
262  {
263  std::stringstream buffer;
264  buffer << "Rigid Body Element #" << this->Id();
265  return buffer.str();
266  }
267 
269  void PrintInfo(std::ostream& rOStream) const override
270  {
271  rOStream << "Rigid Body Element #" << this->Id();
272  }
273 
275  void PrintData(std::ostream& rOStream) const override
276  {
277  GetGeometry().PrintData(rOStream);
278  }
283 
284 protected:
285 
291 
293 
300 
304  virtual void SetProcessInformation(const ProcessInfo& rCurrentProcessInfo);
305 
306 
310  void GetTimeIntegrationParameters(double& rP0,double& rP1,double& rP2,
311  const ProcessInfo& rCurrentProcessInfo) override;
312 
316  SizeType GetDofsSize() const override;
317 
321  void UpdateRigidBodyNodes(const ProcessInfo& rCurrentProcessInfo) override;
322 
323 
334 
335 private:
336 
355  friend class Serializer;
356 
357  virtual void save(Serializer& rSerializer) const override;
358 
359  virtual void load(Serializer& rSerializer) override;
360 
367 
368 
369 }; // Class TranslatoryRigidBodySegregatedVElement
370 
371 } // namespace Kratos.
372 #endif // KRATOS_TRANSLATORY_RIGID_BODY_SEGREGATED_V_ELEMENT_H_INCLUDED defined
Definition: beam_math_utilities.hpp:31
std::size_t SizeType
Definition: element.h:94
std::size_t SizeType
Definition: geometry_data.h:173
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
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
Rigid Body Element for 3D space dimension.
Definition: translatory_rigid_body_element.hpp:46
Rigid Body Segregated V Element for 3D space dimension.
Definition: translatory_rigid_body_segregated_V_element.hpp:46
Node NodeType
Type for nodes.
Definition: translatory_rigid_body_segregated_V_element.hpp:56
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(TranslatoryRigidBodySegregatedVElement)
Counted pointer of TranslatoryRigidBodySegregatedVElement.
StepType mStepVariable
Definition: translatory_rigid_body_segregated_V_element.hpp:292
PointerVectorSet< NodeType, IndexedObject > NodesContainerType
Type for nodes container.
Definition: translatory_rigid_body_segregated_V_element.hpp:58
Quaternion< double > QuaternionType
Type definition for quaternion.
Definition: translatory_rigid_body_segregated_V_element.hpp:54
BeamMathUtils< double > BeamMathUtilsType
Definition: translatory_rigid_body_segregated_V_element.hpp:52
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: translatory_rigid_body_segregated_V_element.hpp:275
TranslatoryRigidBodySegregatedVElement()
Serializer constructor.
Definition: translatory_rigid_body_segregated_V_element.hpp:72
StepType
Definition: translatory_rigid_body_segregated_V_element.hpp:65
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: translatory_rigid_body_segregated_V_element.hpp:269
std::string Info() const override
Turn back information as a string.
Definition: translatory_rigid_body_segregated_V_element.hpp:261
GeometryData::SizeType SizeType
Type for size.
Definition: translatory_rigid_body_segregated_V_element.hpp:60
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
ModelPart::NodesContainerType NodesArrayType
Definition: gid_gauss_point_container.h:42
def load(f)
Definition: ode_solve.py:307