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.
scalar_wall_flux_condition.h
Go to the documentation of this file.
1 // | / |
2 // ' / __| _` | __| _ \ __|
3 // . \ | ( | | ( |\__ `
4 // _|\_\_| \__,_|\__|\___/ ____/
5 // Multi-Physics
6 //
7 // License: BSD License
8 // Kratos default license: kratos/license.txt
9 //
10 // Main authors: Suneth Warnakulasuriya
11 //
12 
13 #if !defined(KRATOS_RANS_SCALAR_WALL_FLUX_CONDITION_H_INCLUDED)
14 #define KRATOS_RANS_SCALAR_WALL_FLUX_CONDITION_H_INCLUDED
15 
16 // System includes
17 
18 // External includes
19 
20 // Project includes
21 #include "includes/condition.h"
22 #include "includes/define.h"
23 
24 // Application includes
25 
26 namespace Kratos
27 {
30 
31 template <unsigned int TDim, unsigned int TNumNodes, class TScalarWallFluxConditionData>
33 {
34 public:
37 
39 
41  using NodeType = Node;
42 
45 
48 
50  using VectorType = Vector;
51 
53  using MatrixType = Matrix;
54 
55  using IndexType = std::size_t;
56 
57  using EquationIdVectorType = std::vector<IndexType>;
58 
59  using DofsVectorType = std::vector<Dof<double>::Pointer>;
60 
61  using ScalarWallFluxConditionDataType = TScalarWallFluxConditionData;
62 
65 
70 
74 
79  IndexType NewId = 0)
80  : Condition(NewId)
81  {
82  }
83 
88  IndexType NewId,
89  const NodesArrayType& ThisNodes)
90  : Condition(NewId, ThisNodes)
91  {
92  }
93 
98  IndexType NewId,
99  GeometryType::Pointer pGeometry)
100  : Condition(NewId, pGeometry)
101  {
102  }
103 
108  IndexType NewId,
109  GeometryType::Pointer pGeometry,
110  PropertiesType::Pointer pProperties)
111  : Condition(NewId, pGeometry, pProperties)
112  {
113  }
114 
119  ScalarWallFluxCondition const& rOther)
120  : Condition(rOther)
121  {
122  }
123 
127  ~ScalarWallFluxCondition() override = default;
128 
132 
145  Condition::Pointer Create(
146  IndexType NewId,
147  NodesArrayType const& ThisNodes,
148  PropertiesType::Pointer pProperties) const override
149  {
150  KRATOS_TRY
151  return Kratos::make_intrusive<CurrentConditionType>(
152  NewId, Condition::GetGeometry().Create(ThisNodes), pProperties);
153  KRATOS_CATCH("");
154  }
155 
163  Condition::Pointer Create(
164  IndexType NewId,
165  GeometryType::Pointer pGeom,
166  PropertiesType::Pointer pProperties) const override
167  {
168  KRATOS_TRY
169  return Kratos::make_intrusive<CurrentConditionType>(NewId, pGeom, pProperties);
170  KRATOS_CATCH("");
171  }
172 
180  Condition::Pointer Clone(
181  IndexType NewId,
182  NodesArrayType const& ThisNodes) const override
183  {
184  KRATOS_TRY
185  return Kratos::make_intrusive<CurrentConditionType>(
187  KRATOS_CATCH("");
188  }
189 
190  void EquationIdVector(
191  EquationIdVectorType& rResult,
192  const ProcessInfo& CurrentProcessInfo) const override;
193 
199  void GetDofList(
200  DofsVectorType& rConditionalDofList,
201  const ProcessInfo& CurrentProcessInfo) const override;
202 
203  void GetValuesVector(
204  Vector& rValues,
205  int Step = 0) const override;
206 
208  Vector& rValues,
209  int Step = 0) const override;
210 
212  Vector& rValues,
213  int Step = 0) const override;
214 
231  MatrixType& rLeftHandSideMatrix,
232  VectorType& rRightHandSideVector,
233  const ProcessInfo& rCurrentProcessInfo) override;
234 
242  VectorType& rRightHandSideVector,
243  const ProcessInfo& rCurrentProcessInfo) override;
244 
254  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
255 
259 
261  std::string Info() const override
262  {
263  std::stringstream buffer;
264  buffer << "ScalarWallFluxCondition #" << Id();
265  return buffer.str();
266  }
268  void PrintInfo(std::ostream& rOStream) const override
269  {
270  rOStream << "SWF" << TScalarWallFluxConditionData::GetName();
271  }
272 
274 
275 private:
278 
279  friend class Serializer;
280 
281  void save(Serializer& rSerializer) const override
282  {
283  KRATOS_TRY
284 
286 
287  KRATOS_CATCH("");
288  }
289  void load(Serializer& rSerializer) override
290  {
291  KRATOS_TRY
292 
294 
295  KRATOS_CATCH("");
296  }
297 
299 
300 }; // Class ScalarWallFluxCondition
301 
305 
306 template <unsigned int TDim, unsigned int TNumNodes, class TScalarWallFluxConditionData>
307 inline std::istream& operator>>(
308  std::istream& rIStream,
310 
312 template <unsigned int TDim, unsigned int TNumNodes, class TScalarWallFluxConditionData>
313 inline std::ostream& operator<<(
314  std::ostream& rOStream,
316 {
317  rThis.PrintInfo(rOStream);
318  rOStream << " : " << std::endl;
319  rThis.PrintData(rOStream);
320  return rOStream;
321 }
322 
324 
325 } // namespace Kratos.
326 
327 #endif // KRATOS_RANS_SCALAR_WALL_FLUX_CONDITION_H_INCLUDED defined
Base class for all Conditions.
Definition: condition.h:59
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: condition.h:1085
std::vector< std::size_t > EquationIdVectorType
Definition: condition.h:98
std::vector< DofType::Pointer > DofsVectorType
Definition: condition.h:100
PropertiesType::Pointer pGetProperties()
returns the pointer to the property of the condition. Does not throw an error, to allow copying of co...
Definition: condition.h:964
Condition(IndexType NewId=0)
Definition: condition.h:123
std::size_t IndexType
Definition: flags.h:74
This defines the geometrical object, base definition of the element and condition entities.
Definition: geometrical_object.h:58
GeometryType & GetGeometry()
Returns the reference of the geometry.
Definition: geometrical_object.h:158
Geometry base class.
Definition: geometry.h:71
IndexType Id() const
Definition: indexed_object.h:107
This class defines the node.
Definition: node.h:65
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
Definition: scalar_wall_flux_condition.h:33
void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: scalar_wall_flux_condition.cpp:116
std::vector< Dof< double >::Pointer > DofsVectorType
Definition: scalar_wall_flux_condition.h:59
Condition::Pointer Create(IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties) const override
Definition: scalar_wall_flux_condition.h:163
ScalarWallFluxCondition(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
Definition: scalar_wall_flux_condition.h:107
std::string Info() const override
Turn back information as a string.
Definition: scalar_wall_flux_condition.h:261
void GetSecondDerivativesVector(Vector &rValues, int Step=0) const override
Definition: scalar_wall_flux_condition.cpp:96
~ScalarWallFluxCondition() override=default
ScalarWallFluxCondition(IndexType NewId, const NodesArrayType &ThisNodes)
Definition: scalar_wall_flux_condition.h:87
void GetDofList(DofsVectorType &rConditionalDofList, const ProcessInfo &CurrentProcessInfo) const override
Definition: scalar_wall_flux_condition.cpp:51
ScalarWallFluxCondition(ScalarWallFluxCondition const &rOther)
Definition: scalar_wall_flux_condition.h:118
void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &CurrentProcessInfo) const override
Definition: scalar_wall_flux_condition.cpp:34
void GetValuesVector(Vector &rValues, int Step=0) const override
Definition: scalar_wall_flux_condition.cpp:68
void GetFirstDerivativesVector(Vector &rValues, int Step=0) const override
Definition: scalar_wall_flux_condition.cpp:76
KRATOS_CLASS_POINTER_DEFINITION(ScalarWallFluxCondition)
std::vector< IndexType > EquationIdVectorType
Definition: scalar_wall_flux_condition.h:57
Matrix MatrixType
Matrix type for local contributions to the linear system.
Definition: scalar_wall_flux_condition.h:53
TScalarWallFluxConditionData ScalarWallFluxConditionDataType
Definition: scalar_wall_flux_condition.h:61
void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: scalar_wall_flux_condition.cpp:133
Condition::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
Definition: scalar_wall_flux_condition.h:145
ScalarWallFluxCondition(IndexType NewId=0)
Definition: scalar_wall_flux_condition.h:78
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: scalar_wall_flux_condition.h:268
ScalarWallFluxCondition(IndexType NewId, GeometryType::Pointer pGeometry)
Definition: scalar_wall_flux_condition.h:97
int Check(const ProcessInfo &rCurrentProcessInfo) const override
Definition: scalar_wall_flux_condition.cpp:176
Condition::Pointer Clone(IndexType NewId, NodesArrayType const &ThisNodes) const override
Definition: scalar_wall_flux_condition.h:180
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
#define KRATOS_SERIALIZE_SAVE_BASE_CLASS(Serializer, BaseType)
Definition: define.h:812
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
Internals::Matrix< double, AMatrix::dynamic, 1 > Vector
Definition: amatrix_interface.h:472
Internals::Matrix< double, AMatrix::dynamic, AMatrix::dynamic > Matrix
Definition: amatrix_interface.h:470
std::istream & operator>>(std::istream &rIStream, LinearMasterSlaveConstraint &rThis)
input stream function
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
def load(f)
Definition: ode_solve.py:307