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_face.h
Go to the documentation of this file.
1 // KRATOS ___ ___ _ ___ __ ___ ___ ___ ___
2 // / __/ _ \| \| \ \ / /__| \_ _| __| __|
3 // | (_| (_) | .` |\ V /___| |) | || _|| _|
4 // \___\___/|_|\_| \_/ |___/___|_| |_| APPLICATION
5 //
6 // License: BSD License
7 // Kratos default license: kratos/license.txt
8 //
9 // Main authors: Jordi Cotela
10 // Riccardo Rossi
11 // Ruben Zorrilla
12 //
13 
14 #if !defined(KRATOS_THERMAL_FACE_H_INCLUDED )
15 #define KRATOS_THERMAL_FACE_H_INCLUDED
16 
17 // System includes
18 #include <string>
19 #include <iostream>
20 
21 // External includes
22 
23 // Project includes
24 #include "includes/define.h"
25 #include "includes/node.h"
26 #include "includes/process_info.h"
27 #include "includes/properties.h"
28 #include "includes/condition.h"
29 #include "geometries/geometry.h"
30 #include "includes/variables.h"
31 
32 namespace Kratos
33 {
34 
37 
39 
43 class KRATOS_API(CONVECTION_DIFFUSION_APPLICATION) ThermalFace: public Condition
44 {
45 public:
48 
51 
57  {
58  double Weight; // Gauss point weight
59  array_1d<double, 3> Normal; // Condition normal
60  Vector N; // Gauss point shape functions values
61 
62  double Emissivity; // Ambient emissivity value
63  double AmbientTemperature; // Ambient temperature value
64  double ConvectionCoefficient; // Ambient convection coefficient
65  Vector UnknownValues; // Previous iteration unknown values
66  Vector FaceHeatFluxValues; // Nodal face heat flux values
67 
68  double inline GaussPointUnknown() const
69  {
70  return InterpolateInGaussPoint(UnknownValues);
71  }
72 
73  double inline GaussPointFaceHeatFlux() const
74  {
75  return InterpolateInGaussPoint(FaceHeatFluxValues);
76  }
77 
78  double inline InterpolateInGaussPoint(const Vector &rNodalValues) const
79  {
80  double gauss_pt_val = 0.0;
81  for (unsigned int i = 0; i < N.size(); ++i) {
82  gauss_pt_val += N[i] * rNodalValues[i];
83  }
84  return gauss_pt_val;
85  }
86  };
87 
90 
92  constexpr static double StefanBoltzmann = 5.67e-8;
93 
97 
99  IndexType NewId,
100  Geometry< Node >::Pointer pGeometry);
101 
102  ThermalFace(
103  IndexType NewId,
104  Geometry< Node >::Pointer pGeometry,
105  Properties::Pointer pProperties);
106 
108  ~ThermalFace() override;
109 
113 
114  Condition::Pointer Create(
115  IndexType NewId,
116  NodesArrayType const& ThisNodes,
117  Properties::Pointer pProperties) const override;
118 
119  Condition::Pointer Create(
120  IndexType NewId,
121  GeometryType::Pointer pGeom,
122  Properties::Pointer pProperties) const override;
123 
124  void CalculateLocalSystem(
125  MatrixType& rLeftHandSideMatrix,
126  VectorType& rRightHandSideVector,
127  const ProcessInfo& rCurrentProcessInfo) override;
128 
129  void CalculateLeftHandSide(
130  MatrixType& rLeftHandSideMatrix,
131  const ProcessInfo& rCurrentProcessInfo) override;
132 
133  void CalculateRightHandSide(
134  VectorType& rRightHandSideVector,
135  const ProcessInfo& rCurrentProcessInfo) override;
136 
137  void EquationIdVector(
138  EquationIdVectorType& rResult,
139  const ProcessInfo& rCurrentProcessInfo) const override;
140 
141  void GetDofList(
142  DofsVectorType& ConditionalDofList,
143  const ProcessInfo& CurrentProcessInfo) const override;
144 
145  GeometryData::IntegrationMethod GetIntegrationMethod() const override;
146 
148  const Variable<array_1d<double, 3 > >& rVariable,
149  std::vector<array_1d<double, 3 > >& rValues,
150  const ProcessInfo& rCurrentProcessInfo) override;
151 
153  const Variable<double>& rVariable,
154  std::vector<double>& rValues,
155  const ProcessInfo& rCurrentProcessInfo) override;
156 
158  const Variable<array_1d<double, 6 > >& rVariable,
159  std::vector<array_1d<double, 6 > >& rValues,
160  const ProcessInfo& rCurrentProcessInfo) override;
161 
163  const Variable<Vector>& rVariable,
164  std::vector<Vector>& rValues,
165  const ProcessInfo& rCurrentProcessInfo) override;
166 
168  const Variable<Matrix>& rVariable,
169  std::vector<Matrix>& rValues,
170  const ProcessInfo& rCurrentProcessInfo) override;
171 
175 
177  std::string Info() const override;
178 
180  void PrintInfo(std::ostream& rOStream) const override;
181 
183  void PrintData(std::ostream& rOStream) const override;
184 
186 
187 protected:
188 
191 
192  // Internal default constructor for serialization
193  ThermalFace();
194 
198 
207  virtual void SetIntegrationWeight(
208  const IndexType IntegrationPointIndex,
209  const typename GeometryType::IntegrationPointsArrayType& rIntegrationPoints,
210  const Vector& rJacobianDeterminantsVector,
211  ConditionDataStruct& rData);
212 
213  void AddIntegrationPointRHSContribution(
214  VectorType& rRightHandSideVector,
215  const ConditionDataStruct &rData);
216 
217  void AddIntegrationPointLHSContribution(
218  MatrixType& rLeftHandSideMatrix,
219  const ConditionDataStruct &rData);
220 
221  void FillConditionDataStructure(
222  const ProcessInfo &rCurrentProcessInfo,
223  ConditionDataStruct &rData);
224 
226 
227 private:
228 
231  friend class Serializer;
232 
233  void save(Serializer& rSerializer) const override;
234 
235  void load(Serializer& rSerializer) override;
236 
240 
242  ThermalFace& operator=(ThermalFace const& rOther);
243 
245  ThermalFace(ThermalFace const& rOther);
246 
248 
249 }; // Class ThermalFace
250 
252 
254 
255 } // namespace Kratos.
256 
257 #endif // KRATOS_THERMAL_FACE_H_INCLUDED defined
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
std::string Info() const override
Turn back information as a string.
Definition: periodic_interface_process.hpp:93
Base class for all Conditions.
Definition: condition.h:59
std::vector< std::size_t > EquationIdVectorType
Definition: condition.h:98
std::vector< DofType::Pointer > DofsVectorType
Definition: condition.h:100
std::size_t IndexType
Definition: flags.h:74
IntegrationMethod
Definition: geometry_data.h:76
Geometry base class.
Definition: geometry.h:71
std::vector< IntegrationPointType > IntegrationPointsArrayType
Definition: geometry.h:161
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
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
A basic Neumann condition for convection-diffusion problems.
Definition: thermal_face.h:44
Condition::MatrixType MatrixType
Definition: thermal_face.h:88
Condition::VectorType VectorType
Definition: thermal_face.h:89
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(ThermalFace)
Pointer definition of ThermalFace.
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
Modeler::Pointer Create(const std::string &ModelerName, Model &rModel, const Parameters ModelParameters)
Checks if the modeler is registered.
Definition: modeler_factory.cpp:30
pybind11::list CalculateOnIntegrationPoints(TObject &dummy, const Variable< TDataType > &rVariable, const ProcessInfo &rProcessInfo)
Definition: add_mesh_to_python.cpp:142
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
def load(f)
Definition: ode_solve.py:307
N
Definition: sensitivityMatrix.py:29
integer i
Definition: TensorModule.f:17
Gauss pt. data structure Auxiliar data structure to pass the Gauss pt. data.
Definition: thermal_face.h:57
double AmbientTemperature
Definition: thermal_face.h:63
double GaussPointUnknown() const
Definition: thermal_face.h:68
double GaussPointFaceHeatFlux() const
Definition: thermal_face.h:73
double Emissivity
Definition: thermal_face.h:62
double ConvectionCoefficient
Definition: thermal_face.h:64
Vector UnknownValues
Definition: thermal_face.h:65
Vector N
Definition: thermal_face.h:60
array_1d< double, 3 > Normal
Definition: thermal_face.h:59
Vector FaceHeatFluxValues
Definition: thermal_face.h:66
double InterpolateInGaussPoint(const Vector &rNodalValues) const
Definition: thermal_face.h:78
double Weight
Definition: thermal_face.h:58