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_penalty_2D_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_PENALTY_2D_CONDITION_H_INCLUDED )
11 #define KRATOS_THERMAL_CONTACT_DOMAIN_PENALTY_2D_CONDITION_H_INCLUDED
12 
13 // System includes
14 
15 // External includes
16 
17 // Project includes
19 
20 namespace Kratos
21 {
36 
37 
38 class KRATOS_API(CONTACT_MECHANICS_APPLICATION) ThermalContactDomainPenalty2DCondition
40 {
41 public:
42 
48  typedef ConstitutiveLawType::Pointer ConstitutiveLawPointerType;
51 
53  typedef Node NodeType;
58 
59 
68 
69 
75 
76 
79 
80  ThermalContactDomainPenalty2DCondition(IndexType NewId, GeometryType::Pointer pGeometry);
81 
82  ThermalContactDomainPenalty2DCondition(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties);
83 
86 
87 
90 
94 
97 
98 
102 
110  Condition::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes, PropertiesType::Pointer pProperties) const override;
111 
112 
113  //************************************************************************************
114  //************************************************************************************
122  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
123 
124  //std::string Info() const;
125 
129 
136 
138  std::string Info() const override
139  {
140  std::stringstream buffer;
141  buffer << "Thermal Contact Domain Penalty 2D Condition #" << Id();
142  return buffer.str();
143 
144  }
146  void PrintInfo(std::ostream& rOStream) const override
147  {
148  rOStream << "Thermal Contact Domain Penalty 2D Condition #" << Id();
149  }
150 
152  void PrintData(std::ostream& rOStream) const override
153  {
154  GetGeometry().PrintData(rOStream);
155  }
156 
161 
162 protected:
171 
175  void SetMasterGeometry() override;
176 
180  void CalculateKinematics(GeneralVariables& rVariables,
181  const ProcessInfo& rCurrentProcessInfo,
182  const unsigned int& rPointNumber) override;
183 
187  void CalcProjections(GeneralVariables& rVariables, const ProcessInfo& rCurrentProcessInfo);
188 
192  double& CalculateIntegrationWeight(double& rIntegrationWeight) override;
193 
194 
198  void CalculateThermalFrictionForce(double &F, GeneralVariables& rVariables, unsigned int& ndi) override;
199 
200 
204  void CalculateThermalConductionForce(double &F, GeneralVariables& rVariables, unsigned int& ndi) override;
205 
206 
210  PointType & CalculateCurrentTangent(PointType &rTangent) override;
211 
221  friend class Serializer;
222 
223  void save(Serializer& rSerializer) const override;
224 
225  void load(Serializer& rSerializer) override;
226 
234 
235 private:
236 
262 
263 }; // Class ThermalContactDomainPenalty2DCondition
264 
272 
273 } // namespace Kratos.
274 #endif // KRATOS_THERMAL_CONTACT_DOMAIN_PENALTY_2D_CONDITION_H_INCLUDED defined
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
Definition: constitutive_law.h:47
Base class for all Elements.
Definition: element.h:60
IntegrationMethod
Definition: geometry_data.h:76
Geometry base class.
Definition: geometry.h:71
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
Definition: thermal_contact_domain_penalty_2D_condition.hpp:40
ContactDomainUtilities::SurfaceVector SurfaceVector
SurfaceVector.
Definition: thermal_contact_domain_penalty_2D_condition.hpp:63
ConstitutiveLaw ConstitutiveLawType
Definition: thermal_contact_domain_penalty_2D_condition.hpp:46
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: thermal_contact_domain_penalty_2D_condition.hpp:146
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(ThermalContactDomainPenalty2DCondition)
Counted pointer of ThermalContactDomainPenalty2DCondition.
std::string Info() const override
Turn back information as a string.
Definition: thermal_contact_domain_penalty_2D_condition.hpp:138
Element::ElementType ElementType
Element Type.
Definition: thermal_contact_domain_penalty_2D_condition.hpp:57
ThermalContactDomainPenalty2DCondition()
Default constructors.
Definition: thermal_contact_domain_penalty_2D_condition.hpp:78
ContactDomainUtilities::SurfaceScalar SurfaceScalar
SurfaceScalar.
Definition: thermal_contact_domain_penalty_2D_condition.hpp:65
GeometryData::IntegrationMethod IntegrationMethod
Type definition for integration methods.
Definition: thermal_contact_domain_penalty_2D_condition.hpp:50
ContactDomainUtilities::PointType PointType
Tensor order 1 definition.
Definition: thermal_contact_domain_penalty_2D_condition.hpp:61
Geometry< NodeType > GeometryType
Geometry Type.
Definition: thermal_contact_domain_penalty_2D_condition.hpp:55
ContactDomainUtilities::BaseLengths BaseLengths
BaseLengths.
Definition: thermal_contact_domain_penalty_2D_condition.hpp:67
ConstitutiveLawType::Pointer ConstitutiveLawPointerType
Pointer type for constitutive laws.
Definition: thermal_contact_domain_penalty_2D_condition.hpp:48
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: thermal_contact_domain_penalty_2D_condition.hpp:152
Node NodeType
NodeType.
Definition: thermal_contact_domain_penalty_2D_condition.hpp:53
std::size_t IndexType
The definition of the index type.
Definition: key_hash.h:35
Modeler::Pointer Create(const std::string &ModelerName, Model &rModel, const Parameters ModelParameters)
Checks if the modeler is registered.
Definition: modeler_factory.cpp:30
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