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_point_link_condition.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: September 2018 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_RIGID_BODY_POINT_LINK_CONDITION_H_INCLUDED )
11 #define KRATOS_RIGID_BODY_POINT_LINK_CONDITION_H_INCLUDED
12 
13 // System includes
14 
15 // External includes
16 
17 // Project includes
18 #include "includes/element.h"
19 #include "includes/condition.h"
21 
22 namespace Kratos
23 {
38 
40 
45 class KRATOS_API(CONTACT_MECHANICS_APPLICATION) RigidBodyPointLinkCondition
46  : public Condition
47 {
48  public:
49 
51 
52  typedef Vector VectorType;
54  typedef Node::Pointer PointPointerType;
59 
62  // Counted pointer of RigidBodyPointLinkCondition
65 
66  protected:
67 
71  KRATOS_DEFINE_LOCAL_FLAG(COMPUTE_RHS_VECTOR);
72  KRATOS_DEFINE_LOCAL_FLAG(COMPUTE_LHS_MATRIX);
76  typedef struct
77  {
78  Flags Options; //calculation options
79 
82 
86 
87  std::vector<SizeType> RigidNodes;
88  std::vector<SizeType> DeformableNodes;
89 
91  std::vector<BoundedMatrix<double,3,3>> RigidSkewSymDistances;
92 
94 
96 
105  {
106  private:
107 
108  //for calculation local system with compacted LHS and RHS
109  MatrixType *mpLeftHandSideMatrix;
110  VectorType *mpRightHandSideVector;
111 
112  public:
113 
114  //calculation flags
116 
120  void SetLeftHandSideMatrix( MatrixType& rLeftHandSideMatrix ) { mpLeftHandSideMatrix = &rLeftHandSideMatrix; };
121  void SetRightHandSideVector( VectorType& rRightHandSideVector ) { mpRightHandSideVector = &rRightHandSideVector; };
122 
126  MatrixType& GetLeftHandSideMatrix() { return *mpLeftHandSideMatrix; };
127  VectorType& GetRightHandSideVector() { return *mpRightHandSideVector; };
128 
129  };
130 
131 
132  public:
133 
134 
137 
139  RigidBodyPointLinkCondition( IndexType NewId, GeometryType::Pointer pGeometry );
140 
141  RigidBodyPointLinkCondition( IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties );
142 
145 
147  virtual ~RigidBodyPointLinkCondition();
148 
152 
153 
157 
165  Condition::Pointer Create(IndexType NewId,
166  NodesArrayType const& ThisNodes,
167  PropertiesType::Pointer pProperties ) const override;
175  Condition::Pointer Clone(IndexType NewId,
176  NodesArrayType const& ThisNodes) const override;
177 
178  //************* STARTING - ENDING METHODS
179 
185  void Initialize(const ProcessInfo& rCurrentProcessInfo) override;
186 
190  void InitializeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo) override;
191 
195  void FinalizeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo) override;
196 
200  void InitializeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
201 
205  void FinalizeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
206 
207  //************* GETTING METHODS
208 
212  void GetDofList(DofsVectorType& rConditionDofList,
213  const ProcessInfo& rCurrentProcessInfo ) const override;
214 
218  void EquationIdVector(EquationIdVectorType& rResult,
219  const ProcessInfo& rCurrentProcessInfo ) const override;
220 
224  void GetValuesVector(Vector& rValues,
225  int Step = 0 ) const override;
226 
230  void GetFirstDerivativesVector(Vector& rValues,
231  int Step = 0 ) const override;
232 
236  void GetSecondDerivativesVector(Vector& rValues,
237  int Step = 0 ) const override;
238 
239 
240  //************* COMPUTING METHODS
241 
250  void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix,
251  VectorType& rRightHandSideVector,
252  const ProcessInfo& rCurrentProcessInfo) override;
253 
260  void CalculateRightHandSide(VectorType& rRightHandSideVector,
261  const ProcessInfo& rCurrentProcessInfo) override;
262 
270  void CalculateSecondDerivativesContributions(MatrixType& rLeftHandSideMatrix,
271  VectorType& rRightHandSideVector,
272  const ProcessInfo& rCurrentProcessInfo) override;
273 
280  void CalculateSecondDerivativesLHS(MatrixType& rLeftHandSideMatrix,
281  const ProcessInfo& rCurrentProcessInfo) override;
282 
283 
290  void CalculateSecondDerivativesRHS(VectorType& rRightHandSideVector,
291  const ProcessInfo& rCurrentProcessInfo) override;
292 
299  void CalculateMassMatrix(MatrixType& rMassMatrix,
300  const ProcessInfo& rCurrentProcessInfo) override;
301 
308  void CalculateDampingMatrix(MatrixType& rDampingMatrix,
309  const ProcessInfo& rCurrentProcessInfo) override;
310 
318  virtual int Check(const ProcessInfo& rCurrentProcessInfo) const override;
319 
320 
334 
335  protected:
341 
342  static const std::array<const VariableData,6> mLinearDofs;
343  static const std::array<const VariableData,3> mAngularDofs;
344 
348 
352 
356  virtual void InitializeSystemMatrices(MatrixType& rLeftHandSideMatrix,
357  VectorType& rRightHandSideVector,
358  const SizeType& rSlaveElementSize,
359  Flags& rCalculationFlags);
363  virtual void InitializeGeneralVariables(GeneralVariables& rVariables,
364  const ProcessInfo& rCurrentProcessInfo);
368  virtual void CalculateConditionSystem(LocalSystemComponents& rLocalSystem,
369  LocalSystemComponents& rLinkedSystem,
370  Element* rSlaveElement,
371  const ProcessInfo& rCurrentProcessInfo);
375  virtual void CalculateAndAddLHS(LocalSystemComponents& rLocalSystem,
376  LocalSystemComponents& rLinkedSystem,
377  GeneralVariables& rVariables);
381  virtual void CalculateAndAddRHS(LocalSystemComponents& rLocalSystem,
382  LocalSystemComponents& rLinkedSystem,
383  GeneralVariables& rVariables);
387  virtual void CalculateAndAddTangent(MatrixType& rLeftHandSideMatrix,
388  MatrixType& rLinkedLeftHandSideMatrix,
389  GeneralVariables& rVariables);
393  virtual void CalculateAndAddForces(VectorType& rRightHandSideVector,
394  VectorType& rLinkedRightHandSideVector,
395  GeneralVariables& rVariables);
396 
400  virtual void CalculateAndAddTangentBeam(MatrixType& rLeftHandSideMatrix,
401  MatrixType& rLinkedLeftHandSideMatrix,
402  GeneralVariables& rVariables);
406  virtual void CalculateAndAddForcesBeam(VectorType& rRightHandSideVector,
407  VectorType& rLinkedRightHandSideVector,
408  GeneralVariables& rVariables);
409 
410 
414  void AssembleLocalLHS(MatrixType& rLeftHandSideMatrix,
415  const MatrixType& rLocalLeftHandSideMatrix,
416  const SizeType& local_index,
417  const SizeType& dofs_size,
418  const SizeType& master_index);
422  void AssembleLocalRHS(VectorType& rRightHandSideVector,
423  const VectorType& rLocalRightHandSideVector,
424  const SizeType& local_index,
425  const SizeType& dofs_size,
426  const SizeType& master_index);
430  virtual SizeType GetDofsSize();
431 
435  void WriteMatrixInRows(std::string MatrixName,
436  const MatrixType& rMatrix);
440 
442 
446 
450 
452 
453  private:
456 
460 
464 
468 
472 
476 
480 
481  friend class Serializer;
482 
483  void save( Serializer& rSerializer ) const override
484  {
486  }
487 
488  void load( Serializer& rSerializer ) override
489  {
491  }
492 
493 
494 }; // class RigidBodyPointLinkCondition.
495 
496 } // namespace Kratos.
497 
498 #endif // KRATOS_RIGID_BODY_POINT_LINK_CONDITION_H_INCLUDED defined
Definition: beam_math_utilities.hpp:31
Base class for all Conditions.
Definition: condition.h:59
std::size_t SizeType
Definition: condition.h:94
Base class for all Elements.
Definition: element.h:60
Definition: flags.h:58
std::size_t SizeType
Definition: geometry_data.h:173
std::vector< std::unique_ptr< Dof< double > >> DofsContainerType
The DoF container type definition.
Definition: node.h:92
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
#define KRATOS_SERIALIZE_SAVE_BASE_CLASS(Serializer, BaseType)
Definition: define.h:812
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
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