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.
thermal_contact_domain_condition.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosMachiningApplication $
3 // Created by: $Author: JMCarbonell $
4 // Last modified by: $Co-Author: $
5 // Date: $Date: October 2013 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_THERMAL_CONTACT_DOMAIN_CONDITION_H_INCLUDED )
11 #define KRATOS_THERMAL_CONTACT_DOMAIN_CONDITION_H_INCLUDED
12 
13 // System includes
14 #include <iostream>
15 #include <iomanip>
16 
17 // External includes
18 #include "boost/smart_ptr.hpp"
19 
20 
21 // Project includes
22 #include "includes/define.h"
23 #include "includes/element.h"
24 #include "includes/serializer.h"
25 #include "includes/condition.h"
27 #include "includes/variables.h"
29 #include "includes/condition.h"
30 
32 
33 namespace Kratos
34 {
49 
50 class KRATOS_API(CONTACT_MECHANICS_APPLICATION) ThermalContactDomainCondition
51  : public Condition
52 {
53 public:
54 
55 
61  typedef ConstitutiveLawType::Pointer ConstitutiveLawPointerType;
64 
66  typedef Node NodeType;
71 
72 
81 
85 
86 protected:
87 
92  typedef struct
93  {
94  Flags Options; //calculation options
95 
96  double ThermalGap; //thermal gap
97 
98 
99  //Geometric variables
101 
102  SurfaceVector CurrentSurface; //Current Normal and Tangent
103  SurfaceVector ReferenceSurface; //Reference Normal and Tangent
104 
105  std::vector<BaseLengths> CurrentBase; //Current Base Lengths variables
106  std::vector<BaseLengths> ReferenceBase; //Reference Base Lengths variables
107 
108 
109  //Axisymmetric
112 
113 
114  //friction:
117 
118 
120 
121 
123  {
124 
125  //The stabilization parameter and penalty parameter
127 
128  //Contact condition conectivities
129  std::vector<unsigned int> nodes;
130  std::vector<unsigned int> order;
131  std::vector<unsigned int> slaves;
132 
133  //Pointer Variables
138 
143  void SetMasterGeometry (GeometryType& rGeometry){ mpMasterGeometry = &rGeometry; }
144  void SetMasterElement (ElementType& rElement){ mpMasterElement = &rElement; }
145  void SetMasterCondition (ConditionType& rCondition){ mpMasterCondition = &rCondition; }
146  void SetMasterNode (NodeType& rNode){ mpMasterNode = &rNode; }
147 
152  GeometryType& GetMasterGeometry() { return (*mpMasterGeometry); }
153  ElementType& GetMasterElement() { return (*mpMasterElement); }
154  ConditionType& GetMasterCondition() { return (*mpMasterCondition); }
155  NodeType& GetMasterNode() { return (*mpMasterNode); }
156 
157 
158  };
159 
160 
161 
162 public:
163 
166 
170 
171 
174 
175  ThermalContactDomainCondition(IndexType NewId, GeometryType::Pointer pGeometry);
176 
177  ThermalContactDomainCondition(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties);
178 
181 
182 
185 
189 
192 
193 
197 
205  Condition::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes, PropertiesType::Pointer pProperties) const override;
206 
207 
214  Condition::Pointer Clone(IndexType NewId, NodesArrayType const& ThisNodes) const override;
215 
216  //************* GETTING METHODS
217 
222  IntegrationMethod GetIntegrationMethod() const override;
223 
227  void GetDofList(DofsVectorType& rConditionalDofList, const ProcessInfo& rCurrentProcessInfo) const override;
228 
232  void EquationIdVector(EquationIdVectorType& rResult, const ProcessInfo& rCurrentProcessInfo) const override;
233 
237  void GetValuesVector(Vector& rValues, int Step = 0) const override;
238 
242  void GetFirstDerivativesVector(Vector& rValues, int Step = 0) const override;
243 
247  void GetSecondDerivativesVector(Vector& rValues, int Step = 0) const override;
248 
249 
250 
251  //on integration points:
263  //SET
267  void SetValuesOnIntegrationPoints(const Variable<double>& rVariable, const std::vector<double>& rValues, const ProcessInfo& rCurrentProcessInfo) override;
271  void SetValuesOnIntegrationPoints(const Variable<Vector>& rVariable, const std::vector<Vector>& rValues, const ProcessInfo& rCurrentProcessInfo) override;
272 
276  void SetValuesOnIntegrationPoints(const Variable<Matrix>& rVariable, const std::vector<Matrix>& rValues, const ProcessInfo& rCurrentProcessInfo) override;
277 
278 
279  //************* STARTING - ENDING METHODS
280 
285  void Initialize(const ProcessInfo& CurrentProcessInfo) override;
286 
287 
291  void InitializeSolutionStep(const ProcessInfo& CurrentProcessInfo) override;
292 
296  void InitializeNonLinearIteration(const ProcessInfo& CurrentProcessInfo) override;
297 
301  void FinalizeSolutionStep(const ProcessInfo& CurrentProcessInfo) override;
302 
303 
304 
305  //************* COMPUTING METHODS
306 
316  void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix, VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo) override;
317 
324  void CalculateRightHandSide(VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo) override;
325 
332  void CalculateLeftHandSide (MatrixType& rLeftHandSideMatrix, const ProcessInfo& rCurrentProcessInfo) override;
333 
334 
335  //on integration points:
339  void CalculateOnIntegrationPoints(const Variable<double>& rVariable, std::vector<double>& rOutput, const ProcessInfo& rCurrentProcessInfo) override;
340 
344  void CalculateOnIntegrationPoints(const Variable<Vector>& rVariable, std::vector<Vector>& rOutput, const ProcessInfo& rCurrentProcessInfo) override;
345 
349  void CalculateOnIntegrationPoints(const Variable<Matrix >& rVariable, std::vector< Matrix >& rOutput, const ProcessInfo& rCurrentProcessInfo) override;
350 
351 
352 
353  //************************************************************************************
354  //************************************************************************************
362  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
363 
364  //std::string Info() const;
365 
369 
376 
378  std::string Info() const override
379  {
380  std::stringstream buffer;
381  buffer << "Thermal Contact Domain Condition #" << Id();
382  return buffer.str();
383 
384  }
386  void PrintInfo(std::ostream& rOStream) const override
387  {
388  rOStream << "Thermal Contact Domain Condition #" << Id();
389  }
390 
392  void PrintData(std::ostream& rOStream) const override
393  {
394  GetGeometry().PrintData(rOStream);
395  }
396 
401 
402 protected:
408 
413 
414 
419 
420 
425 
429 
433  virtual void SetMasterGeometry()
434  {
435  KRATOS_THROW_ERROR( std::invalid_argument, "Calling base class in contact domain", "" );
436  };
437 
441  void CalculateHeatConductivity();
442 
443 
447  virtual void CalculateConditionalSystem(MatrixType& rLeftHandSideMatrix,
448  VectorType& rRightHandSideVector,
449  const ProcessInfo& rCurrentProcessInfo,
450  Flags& rCalculationFlags);
451 
452 
453 
457  virtual void CalculateKinematics(GeneralVariables& rVariables,
458  const ProcessInfo& rCurrentProcessInfo,
459  const unsigned int& rPointNumber)
460  {
461  KRATOS_THROW_ERROR( std::invalid_argument, "Calling base class in contact domain", "" );
462  };
463 
464 
465 
470 
471 
472 
476  virtual void InitializeSystemMatrices(MatrixType& rLeftHandSideMatrix,
477  VectorType& rRightHandSideVector,
478  Flags& rCalculationFlags);
479 
480 
484  virtual void CalculateAndAddLHS(MatrixType& rLeftHandSideMatrix,
485  GeneralVariables& rVariables,
486  double& rIntegrationWeight);
487 
491  virtual void CalculateAndAddRHS(VectorType& rRightHandSideVector,
492  GeneralVariables& rVariables,
493  double& rIntegrationWeight);
494 
498  virtual void CalculateAndAddThermalKm(MatrixType& rLeftHandSideMatrix,
499  GeneralVariables& rVariables,
500  double& rIntegrationWeight);
501 
505  virtual void CalculateAndAddThermalContactForces(VectorType& rRightHandSideVector,
506  GeneralVariables& rVariables,
507  double& rIntegrationWeight);
508 
509 
513  virtual double& CalculateIntegrationWeight(double& rIntegrationWeight);
514 
515 
519  virtual void CalculateThermalFrictionForce(double &F, GeneralVariables& rVariables, unsigned int& ndi)
520  {
521  KRATOS_THROW_ERROR( std::invalid_argument, "Calling base class in contact domain", "" );
522 
523  };
524 
525 
529  virtual void CalculateThermalConductionForce(double &F, GeneralVariables& rVariables, unsigned int& ndi)
530  {
531  KRATOS_THROW_ERROR( std::invalid_argument, "Calling base class in contact domain", "" );
532 
533  };
534 
535 
539  void CalculateRelativeVelocity(GeneralVariables& rVariables, PointType& TangentVelocity, const ProcessInfo& rCurrentProcessInfo);
540 
541 
545  void CalculateRelativeDisplacement(GeneralVariables& rVariables, PointType & TangentDisplacement, const ProcessInfo& rCurrentProcessInfo);
546 
547 
552  {
553  KRATOS_THROW_ERROR( std::invalid_argument, "Calling base class in contact domain", "" );
554 
555  };
556 
566  friend class Serializer;
567 
568  void save(Serializer& rSerializer) const override;
569 
570  void load(Serializer& rSerializer) override;
571 
579 
580 private:
581 
607 
608 }; // Class ThermalContactDomainCondition
609 
617 
618 } // namespace Kratos.
619 #endif // KRATOS_THERMAL_CONTACT_DOMAIN_CONDITION_H_INCLUDED defined
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
Base class for all Conditions.
Definition: condition.h:59
Definition: constitutive_law.h:47
Short class definition.
Definition: contact_domain_utilities.hpp:45
Base class for all Elements.
Definition: element.h:60
Definition: flags.h:58
IntegrationMethod
Definition: geometry_data.h:76
Geometry base class.
Definition: geometry.h:71
This class is a vector which stores global pointers.
Definition: global_pointers_vector.h:61
This class defines the node.
Definition: node.h:65
Point class.
Definition: point.h:59
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
Definition: thermal_contact_domain_condition.hpp:52
GeometryData::IntegrationMethod IntegrationMethod
Type definition for integration methods.
Definition: thermal_contact_domain_condition.hpp:63
std::string Info() const override
Turn back information as a string.
Definition: thermal_contact_domain_condition.hpp:378
ContactDomainUtilities::SurfaceScalar SurfaceScalar
SurfaceScalar.
Definition: thermal_contact_domain_condition.hpp:78
virtual void CalculateThermalFrictionForce(double &F, GeneralVariables &rVariables, unsigned int &ndi)
Definition: thermal_contact_domain_condition.hpp:519
Node NodeType
NodeType.
Definition: thermal_contact_domain_condition.hpp:66
ConstitutiveLawType::Pointer ConstitutiveLawPointerType
Pointer type for constitutive laws.
Definition: thermal_contact_domain_condition.hpp:61
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: thermal_contact_domain_condition.hpp:386
ContactDomainUtilities::BaseLengths BaseLengths
BaseLengths.
Definition: thermal_contact_domain_condition.hpp:80
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(ThermalContactDomainCondition)
Counted pointer of ThermalContactDomainCondition.
ContactVariables mContactVariables
Definition: thermal_contact_domain_condition.hpp:418
ContactDomainUtilities::PointType PointType
Tensor order 1 definition.
Definition: thermal_contact_domain_condition.hpp:74
GlobalPointersVector< Condition > ConditionWeakPtrVectorType
Definition: thermal_contact_domain_condition.hpp:84
Element::ElementType ElementType
Element Type.
Definition: thermal_contact_domain_condition.hpp:70
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: thermal_contact_domain_condition.hpp:392
virtual void CalculateKinematics(GeneralVariables &rVariables, const ProcessInfo &rCurrentProcessInfo, const unsigned int &rPointNumber)
Definition: thermal_contact_domain_condition.hpp:457
virtual void CalculateThermalConductionForce(double &F, GeneralVariables &rVariables, unsigned int &ndi)
Definition: thermal_contact_domain_condition.hpp:529
IntegrationMethod mThisIntegrationMethod
Definition: thermal_contact_domain_condition.hpp:412
ConstitutiveLaw ConstitutiveLawType
Definition: thermal_contact_domain_condition.hpp:59
GlobalPointersVector< Element > ElementWeakPtrVectorType
Definition: thermal_contact_domain_condition.hpp:83
virtual void SetMasterGeometry()
Definition: thermal_contact_domain_condition.hpp:433
Geometry< NodeType > GeometryType
Geometry Type.
Definition: thermal_contact_domain_condition.hpp:68
ContactDomainUtilities mContactUtilities
Definition: thermal_contact_domain_condition.hpp:424
ThermalContactDomainCondition()
Default constructors.
Definition: thermal_contact_domain_condition.hpp:173
ContactDomainUtilities::SurfaceVector SurfaceVector
SurfaceVector.
Definition: thermal_contact_domain_condition.hpp:76
GlobalPointersVector< Node > NodeWeakPtrVectorType
Definition: thermal_contact_domain_condition.hpp:82
virtual PointType & CalculateCurrentTangent(PointType &rTangent)
Definition: thermal_contact_domain_condition.hpp:551
#define KRATOS_THROW_ERROR(ExceptionType, ErrorMessage, MoreInfo)
Definition: define.h:77
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 SetValuesOnIntegrationPoints(TObject &dummy, const Variable< TDataType > &rVariable, const std::vector< TDataType > &values, const ProcessInfo &rCurrentProcessInfo)
Definition: add_mesh_to_python.cpp:185
pybind11::list CalculateOnIntegrationPoints(TObject &dummy, const Variable< TDataType > &rVariable, const ProcessInfo &rProcessInfo)
Definition: add_mesh_to_python.cpp:142
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
F
Definition: hinsberg_optimization.py:144
def load(f)
Definition: ode_solve.py:307
Definition: contact_domain_utilities.hpp:74
Definition: contact_domain_utilities.hpp:91
Definition: contact_domain_utilities.hpp:83
Definition: thermal_contact_domain_condition.hpp:123
GeometryType * mpMasterGeometry
Definition: thermal_contact_domain_condition.hpp:134
void SetMasterGeometry(GeometryType &rGeometry)
Definition: thermal_contact_domain_condition.hpp:143
NodeType * mpMasterNode
Definition: thermal_contact_domain_condition.hpp:137
void SetMasterNode(NodeType &rNode)
Definition: thermal_contact_domain_condition.hpp:146
ElementType * mpMasterElement
Definition: thermal_contact_domain_condition.hpp:135
std::vector< unsigned int > nodes
Definition: thermal_contact_domain_condition.hpp:129
Condition * mpMasterCondition
Definition: thermal_contact_domain_condition.hpp:136
NodeType & GetMasterNode()
Definition: thermal_contact_domain_condition.hpp:155
GeometryType & GetMasterGeometry()
Definition: thermal_contact_domain_condition.hpp:152
std::vector< unsigned int > order
Definition: thermal_contact_domain_condition.hpp:130
std::vector< unsigned int > slaves
Definition: thermal_contact_domain_condition.hpp:131
ElementType & GetMasterElement()
Definition: thermal_contact_domain_condition.hpp:153
void SetMasterCondition(ConditionType &rCondition)
Definition: thermal_contact_domain_condition.hpp:145
void SetMasterElement(ElementType &rElement)
Definition: thermal_contact_domain_condition.hpp:144
ConditionType & GetMasterCondition()
Definition: thermal_contact_domain_condition.hpp:154
double StabilizationFactor
Definition: thermal_contact_domain_condition.hpp:126
Definition: thermal_contact_domain_condition.hpp:93
SurfaceVector ReferenceSurface
Definition: thermal_contact_domain_condition.hpp:103
double ReferenceRadius
Definition: thermal_contact_domain_condition.hpp:111
std::vector< BaseLengths > CurrentBase
Definition: thermal_contact_domain_condition.hpp:105
double RelativeVelocityNorm
Definition: thermal_contact_domain_condition.hpp:115
double CurrentRadius
Definition: thermal_contact_domain_condition.hpp:110
double FrictionForceNorm
Definition: thermal_contact_domain_condition.hpp:116
Flags Options
Definition: thermal_contact_domain_condition.hpp:94
SurfaceVector CurrentSurface
Definition: thermal_contact_domain_condition.hpp:102
Vector ProjectionsVector
Definition: thermal_contact_domain_condition.hpp:100
double ThermalGap
Definition: thermal_contact_domain_condition.hpp:96
std::vector< BaseLengths > ReferenceBase
Definition: thermal_contact_domain_condition.hpp:106