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_element.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosSolidMechanicsApplication $
3 // Created by: $Author: JMCarbonell $
4 // Last modified by: $Co-Author: $
5 // Date: $Date: July 2013 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_THERMAL_ELEMENT_H_INCLUDED )
11 #define KRATOS_THERMAL_ELEMENT_H_INCLUDED
12 
13 // System includes
14 
15 // External includes
16 
17 // Project includes
18 #include "includes/checks.h"
19 #include "includes/element.h"
20 #include "custom_utilities/element_utilities.hpp"
21 
22 namespace Kratos
23 {
38 
39 class KRATOS_API(SOLID_MECHANICS_APPLICATION) ThermalElement
40  : public Element
41 {
42  public:
43 
49  typedef ConstitutiveLawType::Pointer ConstitutiveLawPointerType;
54 
55 
58 
59 
60  protected:
61 
66  KRATOS_DEFINE_LOCAL_FLAG( COMPUTE_RHS_VECTOR );
67  KRATOS_DEFINE_LOCAL_FLAG( COMPUTE_LHS_MATRIX );
68 
74  {
75 
76  //for axisymmetric use only
77  double CurrentRadius;
79 
80  //general thermal variables
81  double DeltaTime;
82  double HeatCapacity;
86 
87  //general mechanical variables
88  double detJ;
91 
92  //variables including all integration points
98 
104  {
105  pDN_De=&rDN_De;
106  };
107 
108  void SetShapeFunctions(const Matrix& rNcontainer)
109  {
110  pNcontainer=&rNcontainer;
111  };
112 
113 
119  {
120  return *pDN_De;
121  };
122 
124  {
125  return *pNcontainer;
126  };
127 
128 
129  };
130 
131 
132  public:
133 
137 
138 
140  ThermalElement(IndexType NewId, GeometryType::Pointer pGeometry);
141  ThermalElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties);
142 
144  ThermalElement(ThermalElement const& rOther);
145 
146 
148  ~ThermalElement() override;
149 
153 
155  ThermalElement& operator=(ThermalElement const& rOther);
156 
157 
161 
169  Element::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes, PropertiesType::Pointer pProperties) const override;
170 
178  //Element::Pointer Clone(IndexType NewId, NodesArrayType const& ThisNodes) const;
179 
180 
181  //************* GETTING METHODS
182 
187  IntegrationMethod GetIntegrationMethod() const override;
188 
192  void GetDofList(DofsVectorType& rElementalDofList, const ProcessInfo& rCurrentProcessInfo) const override;
193 
197  void EquationIdVector(EquationIdVectorType& rResult, const ProcessInfo& rCurrentProcessInfo) const override;
198 
202  void GetValuesVector(Vector& rValues, int Step = 0) const override;
203 
207  void GetFirstDerivativesVector(Vector& rValues, int Step = 0) const override;
208 
212  void GetSecondDerivativesVector(Vector& rValues, int Step = 0) const override;
213 
214 
215  //on integration points:
227  //SET
231  void SetValuesOnIntegrationPoints(const Variable<double>& rVariable, const std::vector<double>& rValues, const ProcessInfo& rCurrentProcessInfo) override;
232 
236  void SetValuesOnIntegrationPoints(const Variable<Vector>& rVariable, const std::vector<Vector>& rValues, const ProcessInfo& rCurrentProcessInfo) override;
237 
241  void SetValuesOnIntegrationPoints(const Variable<Matrix>& rVariable, const std::vector<Matrix>& rValues, const ProcessInfo& rCurrentProcessInfo) override;
242 
243 
244  //************* STARTING - ENDING METHODS
245 
250  void Initialize(const ProcessInfo& rCurrentProcessInfo) override;
251 
252 
256  void InitializeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
257 
258 
262  void FinalizeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
263 
264 
265 
266  //************* COMPUTING METHODS
267 
277  void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix, VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo) override;
278 
285  void CalculateRightHandSide(VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo) override;
286 
293  void CalculateLeftHandSide (MatrixType& rLeftHandSideMatrix, const ProcessInfo& rCurrentProcessInfo) override;
294 
301  void CalculateMassMatrix(MatrixType& rMassMatrix, const ProcessInfo& rCurrentProcessInfo) override;
302 
309  void CalculateDampingMatrix(MatrixType& rDampingMatrix, const ProcessInfo& rCurrentProcessInfo) override;
310 
311 
312  //on integration points:
316  void CalculateOnIntegrationPoints(const Variable<double>& rVariable, std::vector<double>& rOutput, const ProcessInfo& rCurrentProcessInfo) override;
317 
321  void CalculateOnIntegrationPoints(const Variable<Vector>& rVariable, std::vector<Vector>& rOutput, const ProcessInfo& rCurrentProcessInfo) override;
322 
326  void CalculateOnIntegrationPoints(const Variable<Matrix >& rVariable, std::vector< Matrix >& rOutput, const ProcessInfo& rCurrentProcessInfo) override;
327 
328 
329  //************************************************************************************
330  //************************************************************************************
338  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
339 
340  //std::string Info() const;
341 
345 
352 
354  // virtual String Info() const;
355 
357  // virtual void PrintInfo(std::ostream& rOStream) const;
358 
360  // virtual void PrintData(std::ostream& rOStream) const;
365 
366  protected:
372 
377 
382  {
383  }
384 
390  virtual void CalculateElementalSystem(MatrixType& rLeftHandSideMatrix,
391  VectorType& rRightHandSideVector,
392  const ProcessInfo& rCurrentProcessInfo,
393  Flags & rCalculationFlags);
397 
398 
403  virtual void CalculateAndAddLHS(MatrixType& rLeftHandSideMatrix,
404  GeneralVariables& rVariables,
405  double& rIntegrationWeight);
406 
411  virtual void CalculateAndAddRHS(VectorType& rRightHandSideVector,
412  GeneralVariables& rVariables,
413  double& rHeatSource,
414  double& rIntegrationWeight);
415 
416 
417 
421  virtual void CalculateAndAddKthermal(MatrixType& rK,
422  GeneralVariables & rVariables,
423  double& rIntegrationWeight
424  );
425 
429  virtual void CalculateAndAddHthermal(MatrixType& rK,
430  GeneralVariables & rVariables,
431  double& rIntegrationWeight
432  );
433 
434 
438  virtual void CalculateAndAddMthermal(MatrixType& rM,
439  GeneralVariables & rVariables,
440  double& rIntegrationWeight
441  );
442 
443 
447  virtual void CalculateAndAddExternalForces(GeneralVariables& rVariables,
448  VectorType& rRightHandSideVector,
449  double& rHeatSource,
450  double& rIntegrationWeight
451  );
452 
453 
457  virtual void CalculateAndAddThermalForces(GeneralVariables & rVariables,
458  VectorType& rRightHandSideVector,
459  double& rIntegrationWeight
460  );
461 
462 
466  virtual void InitializeSystemMatrices(MatrixType& rLeftHandSideMatrix,
467  VectorType& rRightHandSideVector,
468  Flags& rCalculationFlags);
469 
474 
478  virtual void CalculateKinematics(GeneralVariables& rVariables,
479  const double& rPointNumber);
480 
481 
485  void CalculateThermalProperties(GeneralVariables& rVariables);
486 
487 
488 
492  virtual double& CalculateIntegrationWeight(double & rIntegrationWeight);
493 
494 
498  virtual void InitializeGeneralVariables(GeneralVariables & rVariables, const ProcessInfo& rCurrentProcessInfo);
499 
500 
501 
512 
513  private:
514 
520 
521 
525 
526 
530 
531 
536 
540  friend class Serializer;
541 
542  void save(Serializer& rSerializer) const override;
543 
544  void load(Serializer& rSerializer) override;
545 
552 
553 }; // Class ThermalElement
554 
562 
563 } // namespace Kratos.
564 #endif // KRATOS_THERMAL_ELEMENT_H_INCLUDED defined
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
Definition: constitutive_law.h:47
StressMeasure
Definition: constitutive_law.h:69
Base class for all Elements.
Definition: element.h:60
Definition: flags.h:58
IntegrationMethod
Definition: geometry_data.h:76
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_element.hpp:41
ConstitutiveLawType::StressMeasure StressMeasureType
StressMeasure from constitutive laws.
Definition: thermal_element.hpp:51
ThermalElement()
Definition: thermal_element.hpp:381
IntegrationMethod mThisIntegrationMethod
Definition: thermal_element.hpp:376
ConstitutiveLawType::Pointer ConstitutiveLawPointerType
Pointer type for constitutive laws.
Definition: thermal_element.hpp:49
KRATOS_DEFINE_LOCAL_FLAG(COMPUTE_LHS_MATRIX)
GeometryData::IntegrationMethod IntegrationMethod
Type definition for integration methods.
Definition: thermal_element.hpp:53
KRATOS_DEFINE_LOCAL_FLAG(COMPUTE_RHS_VECTOR)
ConstitutiveLaw ConstitutiveLawType
Definition: thermal_element.hpp:47
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(ThermalElement)
Counted pointer of ThermalElement.
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
def load(f)
Definition: ode_solve.py:307
Definition: thermal_element.hpp:74
double PlasticDissipation
Definition: thermal_element.hpp:84
double detJ
Definition: thermal_element.hpp:88
Vector N
Definition: thermal_element.hpp:89
double HeatConductivity
Definition: thermal_element.hpp:83
void SetShapeFunctionsGradients(const GeometryType::ShapeFunctionsGradientsType &rDN_De)
Definition: thermal_element.hpp:103
const Matrix * pNcontainer
Definition: thermal_element.hpp:96
const GeometryType::ShapeFunctionsGradientsType * pDN_De
Definition: thermal_element.hpp:93
GeometryType::JacobiansType j
Definition: thermal_element.hpp:95
Matrix DeltaPosition
Definition: thermal_element.hpp:97
const GeometryType::ShapeFunctionsGradientsType & GetShapeFunctionsGradients()
Definition: thermal_element.hpp:118
double CurrentRadius
Definition: thermal_element.hpp:77
const Matrix & GetShapeFunctions()
Definition: thermal_element.hpp:123
double ReferenceRadius
Definition: thermal_element.hpp:78
Matrix DN_DX
Definition: thermal_element.hpp:90
double DeltaPlasticDissipation
Definition: thermal_element.hpp:85
double DeltaTime
Definition: thermal_element.hpp:81
GeometryType::JacobiansType J
Definition: thermal_element.hpp:94
void SetShapeFunctions(const Matrix &rNcontainer)
Definition: thermal_element.hpp:108
double HeatCapacity
Definition: thermal_element.hpp:82