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.
point_rigid_contact_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: July 2016 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_POINT_RIGID_CONTACT_CONDITION_H_INCLUDED )
11 #define KRATOS_POINT_RIGID_CONTACT_CONDITION_H_INCLUDED
12 
13 // System includes
14 
15 // External includes
16 #include "boost/smart_ptr.hpp"
17 
18 // Project includes
19 #include "includes/define.h"
20 #include "includes/serializer.h"
21 #include "includes/condition.h"
23 #include "includes/variables.h"
24 
27 
28 namespace Kratos
29 {
44 
46 
51 class KRATOS_API(CONTACT_MECHANICS_APPLICATION) PointRigidContactCondition
52  : public Condition
53 {
54 public:
55 
57 
59  //typedef BoundedVector<double, 3> PointType;
65  // Counted pointer of PointRigidContactCondition
68 
69 protected:
70 
75  KRATOS_DEFINE_LOCAL_FLAG( COMPUTE_RHS_VECTOR );
76  KRATOS_DEFINE_LOCAL_FLAG( COMPUTE_LHS_MATRIX );
77  KRATOS_DEFINE_LOCAL_FLAG( COMPUTE_RHS_VECTOR_WITH_COMPONENTS );
78  KRATOS_DEFINE_LOCAL_FLAG( COMPUTE_LHS_MATRIX_WITH_COMPONENTS );
79 
80 
85  typedef struct
86  {
87  PointType Normal; //normal direction
88  PointType Tangent; //tangent direction
89 
90  } SurfaceVector;
91 
92  typedef struct
93  {
94  double Normal; //normal component
95  double Tangent; //tangent component
96 
97  } SurfaceScalar;
98 
100  {
101  Flags Options; //calculation options
102 
103  //Geometrical gaps:
104  SurfaceScalar Gap; //normal and tangential gap
105  double ContributoryFactor; //distance to neighbour points
106 
107  //Friction:
108  bool Slip; //slip = true / stick = false
109  PointType RelativeDisplacement; //relative point displacement
110  double FrictionCoefficient; //total friction coeffitient mu
111  double DeltaTime; //time step
112 
113  SurfaceScalar Penalty; //Penalty Parameter normal and tangent
114 
115  //Geometric variables
116  SurfaceVector Surface; //normal and tangent vector to the surface
117 
118  // Elasto-plastic constitutive matrix (axisym and PS) (Normal and Tangent stiffness)
120 
121  //Contact stress
123 
124  //for axisymmetric use only
127 
128  //for rigid body
129  Matrix SkewSymDistance; //compute the skewsymmmetric tensor of the distance
130 
131  void Initialize()
132  {
133  Gap.Normal = 0;
134  Gap.Tangent = 0;
135 
136  ContributoryFactor = 0;
137 
138  RelativeDisplacement.resize(3);
139  noalias(RelativeDisplacement) = ZeroVector(3);
140 
141  Slip = false;
142  DeltaTime = 0;
143  FrictionCoefficient = 0;
144  Penalty.Normal = 0;
145  Penalty.Tangent = 0;
146 
147  Surface.Normal.resize(3);
148  noalias(Surface.Normal) = ZeroVector(3);
149  Surface.Tangent.resize(3);
150  noalias(Surface.Tangent) = ZeroVector(3);
151 
152  SkewSymDistance.resize(3,3,false);
153  noalias(SkewSymDistance) = ZeroMatrix(3,3);
154 
155  TangentMatrix.Normal = 0;
156  TangentMatrix.Tangent = 0;
157 
158  CurrentRadius = 0;
159  ReferenceRadius = 0;
160  }
161 
162  };
163 
164 
173  {
174  private:
175 
176  //for calculation local system with compacted LHS and RHS
177  MatrixType *mpLeftHandSideMatrix;
178  VectorType *mpRightHandSideVector;
179 
180  //for calculation local system with LHS and RHS components
181  std::vector<MatrixType> *mpLeftHandSideMatrices;
182  std::vector<VectorType> *mpRightHandSideVectors;
183 
184  //LHS variable components
185  const std::vector< Variable< MatrixType > > *mpLeftHandSideVariables;
186 
187  //RHS variable components
188  const std::vector< Variable< VectorType > > *mpRightHandSideVariables;
189 
190 
191  public:
192 
193  //calculation flags
195 
199  void SetLeftHandSideMatrix( MatrixType& rLeftHandSideMatrix ) { mpLeftHandSideMatrix = &rLeftHandSideMatrix; };
200  void SetLeftHandSideMatrices( std::vector<MatrixType>& rLeftHandSideMatrices ) { mpLeftHandSideMatrices = &rLeftHandSideMatrices; };
201  void SetLeftHandSideVariables(const std::vector< Variable< MatrixType > >& rLeftHandSideVariables ) { mpLeftHandSideVariables = &rLeftHandSideVariables; };
202 
203  void SetRightHandSideVector( VectorType& rRightHandSideVector ) { mpRightHandSideVector = &rRightHandSideVector; };
204  void SetRightHandSideVectors( std::vector<VectorType>& rRightHandSideVectors ) { mpRightHandSideVectors = &rRightHandSideVectors; };
205  void SetRightHandSideVariables(const std::vector< Variable< VectorType > >& rRightHandSideVariables ) { mpRightHandSideVariables = &rRightHandSideVariables; };
206 
207 
211  MatrixType& GetLeftHandSideMatrix() { return *mpLeftHandSideMatrix; };
212  std::vector<MatrixType>& GetLeftHandSideMatrices() { return *mpLeftHandSideMatrices; };
213  const std::vector< Variable< MatrixType > >& GetLeftHandSideVariables() { return *mpLeftHandSideVariables; };
214 
215  VectorType& GetRightHandSideVector() { return *mpRightHandSideVector; };
216  std::vector<VectorType>& GetRightHandSideVectors() { return *mpRightHandSideVectors; };
217  const std::vector< Variable< VectorType > >& GetRightHandSideVariables() { return *mpRightHandSideVariables; };
218 
219  };
220 
221 
222 public:
223 
224 
227 
230 
232  PointRigidContactCondition( IndexType NewId, GeometryType::Pointer pGeometry );
233 
234  PointRigidContactCondition( IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties );
235 
236 
237  PointRigidContactCondition( IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties, SpatialBoundingBox::Pointer pRigidWall );
238 
241 
243  virtual ~PointRigidContactCondition();
244 
248 
249 
253 
261  Condition::Pointer Create(IndexType NewId,
262  NodesArrayType const& ThisNodes,
263  PropertiesType::Pointer pProperties ) const override;
264 
265 
273  Condition::Pointer Clone(IndexType NewId,
274  NodesArrayType const& ThisNodes) const override;
275 
276 
277  //************* STARTING - ENDING METHODS
278 
282  void Initialize(const ProcessInfo& rCurrentProcessInfo) override;
283 
287  void InitializeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
288 
292  void InitializeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo) override;
293 
297  void FinalizeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo) override;
298 
302  void FinalizeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
303 
304 
305  //************* GETTING METHODS
306 
310  void GetDofList(DofsVectorType& rConditionDofList,
311  const ProcessInfo& rCurrentProcessInfo ) const override;
312 
316  void EquationIdVector(EquationIdVectorType& rResult,
317  const ProcessInfo& rCurrentProcessInfo ) const override;
318 
322  void GetValuesVector(Vector& rValues,
323  int Step = 0 ) const override;
324 
328  void GetFirstDerivativesVector(Vector& rValues,
329  int Step = 0 ) const override;
330 
334  void GetSecondDerivativesVector(Vector& rValues,
335  int Step = 0 ) const override;
336 
337 
338  //************* COMPUTING METHODS
339 
348  void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix,
349  VectorType& rRightHandSideVector,
350  const ProcessInfo& rCurrentProcessInfo ) override;
351 
352 
362  void CalculateLocalSystem(std::vector< MatrixType >& rLeftHandSideMatrices,
363  const std::vector< Variable< MatrixType > >& rLHSVariables,
364  std::vector< VectorType >& rRightHandSideVectors,
365  const std::vector< Variable< VectorType > >& rRHSVariables,
366  const ProcessInfo& rCurrentProcessInfo);
367 
374  void CalculateRightHandSide(VectorType& rRightHandSideVector,
375  const ProcessInfo& rCurrentProcessInfo ) override;
376 
377 
385  void CalculateRightHandSide(std::vector< VectorType >& rRightHandSideVectors,
386  const std::vector< Variable< VectorType > >& rRHSVariables,
387  const ProcessInfo& rCurrentProcessInfo);
388 
395  void CalculateMassMatrix(
396  MatrixType& rMassMatrix,
397  const ProcessInfo& rCurrentProcessInfo ) override;
398 
405  void CalculateDampingMatrix(
406  MatrixType& rDampingMatrix,
407  const ProcessInfo& rCurrentProcessInfo ) override;
408 
409 
419  void AddExplicitContribution(const VectorType& rRHSVector,
420  const Variable<VectorType>& rRHSVariable,
421  const Variable<array_1d<double,3> >& rDestinationVariable,
422  const ProcessInfo& rCurrentProcessInfo) override;
423 
424  //************************************************************************************
425  //************************************************************************************
433  int Check( const ProcessInfo& rCurrentProcessInfo ) const override;
434 
435 
440  void SetRigidWall(SpatialBoundingBox::Pointer pRigidWall);
441 
442 
456 
457 protected:
463 
468 
469 
473  SpatialBoundingBox::Pointer mpRigidWall;
474 
475 
479  FrictionLaw::Pointer mpFrictionLaw;
480 
481 
486 
493 
497  void ClearNodalForces ();
498 
499 
503  virtual void InitializeSystemMatrices(MatrixType& rLeftHandSideMatrix,
504  VectorType& rRightHandSideVector,
505  Flags& rCalculationFlags);
506 
510  virtual void InitializeConditionVariables(ConditionVariables& rVariables,
511  const ProcessInfo& rCurrentProcessInfo);
512 
516  virtual void CalculateKinematics(ConditionVariables& rVariables,
517  const ProcessInfo& rCurrentProcessInfo,
518  const double& rPointNumber);
519 
523  virtual double& CalculateIntegrationWeight(double& rIntegrationWeight);
524 
525 
529  virtual void CalculateConditionSystem(LocalSystemComponents& rLocalSystem,
530  const ProcessInfo& rCurrentProcessInfo);
531 
532 
536  virtual void CalculateAndAddLHS(LocalSystemComponents& rLocalSystem,
537  ConditionVariables& rVariables,
538  double& rIntegrationWeight);
539 
543  virtual void CalculateAndAddRHS(LocalSystemComponents& rLocalSystem,
544  ConditionVariables& rVariables,
545  double& rIntegrationWeight);
546 
547 
551  virtual void CalculateAndAddKuug(MatrixType& rLeftHandSideMatrix,
552  ConditionVariables& rVariables,
553  double& rIntegrationWeight);
554 
555 
559  virtual void CalculateAndAddContactForces(Vector& rRightHandSideVector,
560  ConditionVariables& rVariables,
561  double& rIntegrationWeight );
562 
566 
567 
571 
572 
576 
577 
579 
580 private:
583 
584 
588 
592 
593 
597 
598 
602 
603 
607 
611 
612  friend class Serializer;
613 
614  void save( Serializer& rSerializer ) const override
615  {
617  //rSerializer.save("mpRigidWall",mpRigidWall); //rebuild in contact search and in restart
618  rSerializer.save("mpFrictionLaw",mpFrictionLaw);
619  rSerializer.save("mContactStressVector",mContactStressVector);
620  }
621 
622  void load( Serializer& rSerializer ) override
623  {
625  //rSerializer.load("mpRigidWall",mpRigidWall);
626  rSerializer.load("mpFrictionLaw",mpFrictionLaw);
627  rSerializer.load("mContactStressVector",mContactStressVector);
628  }
629 
631 
632 }; // class PointRigidContactCondition.
633 
634 } // namespace Kratos.
635 
636 #endif // KRATOS_POINT_RIGID_CONTACT_CONDITION_H_INCLUDED defined
Base class for all Conditions.
Definition: condition.h:59
Definition: flags.h:58
IntegrationMethod
Definition: geometry_data.h:76
This class is a vector which stores global pointers.
Definition: global_pointers_vector.h:61
void resize(std::size_t NewSize1, std::size_t NewSize2, bool preserve=0)
Definition: amatrix_interface.h:224
Point Rigid Contact Condition for 3D and 2D geometries. (base class)
Definition: point_rigid_contact_condition.hpp:53
GlobalPointersVector< Condition > ConditionWeakPtrVectorType
Definition: point_rigid_contact_condition.hpp:63
KRATOS_DEFINE_LOCAL_FLAG(COMPUTE_LHS_MATRIX)
IntegrationMethod mThisIntegrationMethod
Definition: point_rigid_contact_condition.hpp:467
void SetRigidWall(SpatialBoundingBox::Pointer pRigidWall)
FrictionLaw::Pointer mpFrictionLaw
Definition: point_rigid_contact_condition.hpp:479
GlobalPointersVector< Element > ElementWeakPtrVectorType
Definition: point_rigid_contact_condition.hpp:62
KRATOS_DEFINE_LOCAL_FLAG(COMPUTE_RHS_VECTOR_WITH_COMPONENTS)
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(PointRigidContactCondition)
PointRigidContactCondition()
Serialization constructor.
Definition: point_rigid_contact_condition.hpp:229
GlobalPointersVector< Node > NodeWeakPtrVectorType
Definition: point_rigid_contact_condition.hpp:61
KRATOS_DEFINE_LOCAL_FLAG(COMPUTE_RHS_VECTOR)
Vector mContactStressVector
Definition: point_rigid_contact_condition.hpp:485
KRATOS_DEFINE_LOCAL_FLAG(COMPUTE_LHS_MATRIX_WITH_COMPONENTS)
array_1d< double, 3 > PointType
Tensor order 1 definition.
Definition: point_rigid_contact_condition.hpp:60
SpatialBoundingBox::Pointer mpRigidWall
Definition: point_rigid_contact_condition.hpp:473
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
void load(std::string const &rTag, TDataType &rObject)
Definition: serializer.h:207
void save(std::string const &rTag, std::array< TDataType, TDataSize > const &rObject)
Definition: serializer.h:545
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
BOOST_UBLAS_INLINE void resize(size_type array_size, bool preserve=true)
Definition: array_1d.h:242
#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
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
def load(f)
Definition: ode_solve.py:307
Definition: point_rigid_contact_condition.hpp:100
SurfaceVector Surface
Definition: point_rigid_contact_condition.hpp:116
Vector ContactStressVector
Definition: point_rigid_contact_condition.hpp:122
double ReferenceRadius
Definition: point_rigid_contact_condition.hpp:126
Flags Options
Definition: point_rigid_contact_condition.hpp:101
double ContributoryFactor
Definition: point_rigid_contact_condition.hpp:105
Matrix SkewSymDistance
Definition: point_rigid_contact_condition.hpp:129
double DeltaTime
Definition: point_rigid_contact_condition.hpp:111
void Initialize()
Definition: point_rigid_contact_condition.hpp:131
double CurrentRadius
Definition: point_rigid_contact_condition.hpp:125
SurfaceScalar Gap
Definition: point_rigid_contact_condition.hpp:104
SurfaceScalar Penalty
Definition: point_rigid_contact_condition.hpp:113
bool Slip
Definition: point_rigid_contact_condition.hpp:108
PointType RelativeDisplacement
Definition: point_rigid_contact_condition.hpp:109
SurfaceScalar TangentMatrix
Definition: point_rigid_contact_condition.hpp:119
double FrictionCoefficient
Definition: point_rigid_contact_condition.hpp:110
Definition: point_rigid_contact_condition.hpp:173
Flags CalculationFlags
Definition: point_rigid_contact_condition.hpp:194
const std::vector< Variable< VectorType > > & GetRightHandSideVariables()
Definition: point_rigid_contact_condition.hpp:217
void SetRightHandSideVectors(std::vector< VectorType > &rRightHandSideVectors)
Definition: point_rigid_contact_condition.hpp:204
void SetLeftHandSideVariables(const std::vector< Variable< MatrixType > > &rLeftHandSideVariables)
Definition: point_rigid_contact_condition.hpp:201
const std::vector< Variable< MatrixType > > & GetLeftHandSideVariables()
Definition: point_rigid_contact_condition.hpp:213
void SetRightHandSideVector(VectorType &rRightHandSideVector)
Definition: point_rigid_contact_condition.hpp:203
void SetRightHandSideVariables(const std::vector< Variable< VectorType > > &rRightHandSideVariables)
Definition: point_rigid_contact_condition.hpp:205
VectorType & GetRightHandSideVector()
Definition: point_rigid_contact_condition.hpp:215
std::vector< MatrixType > & GetLeftHandSideMatrices()
Definition: point_rigid_contact_condition.hpp:212
MatrixType & GetLeftHandSideMatrix()
Definition: point_rigid_contact_condition.hpp:211
void SetLeftHandSideMatrices(std::vector< MatrixType > &rLeftHandSideMatrices)
Definition: point_rigid_contact_condition.hpp:200
std::vector< VectorType > & GetRightHandSideVectors()
Definition: point_rigid_contact_condition.hpp:216
void SetLeftHandSideMatrix(MatrixType &rLeftHandSideMatrix)
Definition: point_rigid_contact_condition.hpp:199
Definition: point_rigid_contact_condition.hpp:93
double Tangent
Definition: point_rigid_contact_condition.hpp:95
double Normal
Definition: point_rigid_contact_condition.hpp:94
Definition: point_rigid_contact_condition.hpp:86
PointType Tangent
Definition: point_rigid_contact_condition.hpp:88
PointType Normal
Definition: point_rigid_contact_condition.hpp:87