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.
two_fluid_navier_stokes_wall_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: Ruben Zorrilla
11 //
12 
13 #pragma once
14 
15 // System includes
16 #include <string>
17 #include <iostream>
18 
19 
20 // External includes
21 
22 
23 // Project includes
24 
25 
26 // Application includes
28 
29 namespace Kratos
30 {
33 
36 
37 
41 
42 
46 
47 
51 
52 
56 
58 
64 template< unsigned int TDim, unsigned int TNumNodes, class... TWallModel>
65 class KRATOS_API(FLUID_DYNAMICS_APPLICATION) TwoFluidNavierStokesWallCondition : public NavierStokesWallCondition<TDim, TNumNodes, TWallModel...>
66 {
67 public:
70 
73 
74  typedef NavierStokesWallCondition<TDim, TNumNodes, TWallModel...> BaseType;
75 
76  using BaseType::LocalSize;
77 
78  typedef typename BaseType::ConditionDataStruct ConditionDataStruct;
79 
80  typedef Node NodeType;
81 
83 
85 
87 
88  typedef Vector VectorType;
89 
90  typedef Matrix MatrixType;
91 
92  typedef std::size_t IndexType;
93 
94  typedef std::vector<std::size_t> EquationIdVectorType;
95 
96  typedef std::vector< Dof<double>::Pointer > DofsVectorType;
97 
101 
103 
109  : NavierStokesWallCondition<TDim, TNumNodes, TWallModel...>(NewId)
110  {
111  }
112 
114 
119  IndexType NewId,
120  const NodesArrayType& ThisNodes)
121  : NavierStokesWallCondition<TDim, TNumNodes, TWallModel...>(NewId, ThisNodes)
122  {
123  }
124 
126 
131  IndexType NewId,
132  GeometryType::Pointer pGeometry)
133  : NavierStokesWallCondition<TDim, TNumNodes, TWallModel...>(NewId, pGeometry)
134  {
135  }
136 
138 
144  IndexType NewId,
145  GeometryType::Pointer pGeometry,
146  PropertiesType::Pointer pProperties)
147  : NavierStokesWallCondition<TDim, TNumNodes, TWallModel...>(NewId, pGeometry, pProperties)
148  {
149  }
150 
153  NavierStokesWallCondition<TDim, TNumNodes, TWallModel...>(rOther)
154  {
155  }
156 
159 
160 
164 
167  {
168  Condition::operator=(rOther);
169  return *this;
170  }
171 
175 
177 
182  Condition::Pointer Create(
183  IndexType NewId,
184  NodesArrayType const& ThisNodes,
185  PropertiesType::Pointer pProperties) const override
186  {
187  return Kratos::make_intrusive<TwoFluidNavierStokesWallCondition>(NewId, BaseType::GetGeometry().Create(ThisNodes), pProperties);
188  }
189 
191 
196  Condition::Pointer Create(
197  IndexType NewId,
198  GeometryType::Pointer pGeom,
199  PropertiesType::Pointer pProperties) const override
200  {
201  return Kratos::make_intrusive< TwoFluidNavierStokesWallCondition >(NewId, pGeom, pProperties);
202  }
203 
210  Condition::Pointer Clone(
211  IndexType NewId,
212  NodesArrayType const& rThisNodes) const override
213  {
214  Condition::Pointer pNewCondition = Create(NewId, BaseType::GetGeometry().Create( rThisNodes ), BaseType::pGetProperties() );
215 
216  pNewCondition->SetData(this->GetData());
217  pNewCondition->SetFlags(this->GetFlags());
218 
219  return pNewCondition;
220  }
221 
228  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
229 
233 
234 
238 
239 
243 
245  std::string Info() const override
246  {
247  std::stringstream buffer;
248  buffer << "TwoFluidNavierStokesWallCondition" << TDim << "D";
249  return buffer.str();
250  }
251 
253  void PrintInfo(std::ostream& rOStream) const override
254  {
255  rOStream << "TwoFluidNavierStokesWallCondition";
256  }
257 
259  void PrintData(std::ostream& rOStream) const override {}
260 
264 
265 
267 protected:
270 
271 
275 
276 
280 
281 
285 
286  // /**
287  // * @brief Computes the right-hand side of the Navier slip contribution as e.g. described in BEHR2004
288  // * The (Navier) slip length is read as a nodal variable.
289  // * If a smaller value is set, tangential velocities lead to a higher tangential traction.
290  // * Though only tangential velocities should appear, a tangetial projection is added.
291  // * (Reference BEHR2004: https://onlinelibrary.wiley.com/doi/abs/10.1002/fld.663)
292  // * @param rRightHandSideVector reference to the RHS vector
293  // * @param rDataStruct reference to a struct to hand over data
294  // */
295  // void ComputeGaussPointNavierSlipRHSContribution(
296  // array_1d<double,LocalSize>& rRightHandSideVector,
297  // const ConditionDataStruct& rDataStruct) override;
298 
299  // /**
300  // * @brief Computes the left-hand side of the Navier slip contribution as e.g. described in BEHR2004
301  // * The (Navier) slip length is read as a nodal variable.
302  // * If a smaller value is set, tangential velocities lead to a higher tangential traction.
303  // * Though only tangential velocities should appear, a tangetial projection is added.
304  // * (Reference BEHR2004: https://onlinelibrary.wiley.com/doi/abs/10.1002/fld.663)
305  // * @param rLeftHandSideMatrix reference to the LHS matrix
306  // * @param rDataStruct reference to a struct to hand over data
307  // */
308  // void ComputeGaussPointNavierSlipLHSContribution(
309  // BoundedMatrix<double,LocalSize,LocalSize>& rLeftHandSideMatrix,
310  // const ConditionDataStruct& rDataStruct ) override;
311 
315 
316 
320 
321 
325 
326 
328 private:
331 
332 
336 
337 
341 
342  friend class Serializer;
343 
344  void save(Serializer& rSerializer) const override
345  {
347  }
348 
349  void load(Serializer& rSerializer) override
350  {
351  KRATOS_SERIALIZE_LOAD_BASE_CLASS(rSerializer, BaseType);
352  }
353 
357 
358 
362 
363 
367 
368 
372 
373 
377 
378 
380 }; // Class TwoFluidNavierStokesWallCondition
381 
385 
386 
390 
392 template< unsigned int TDim, unsigned int TNumNodes, class TWallModel >
393 inline std::istream& operator >> (std::istream& rIStream, TwoFluidNavierStokesWallCondition<TDim,TNumNodes,TWallModel>& rThis)
394 {
395  return rIStream;
396 }
397 
399 template< unsigned int TDim, unsigned int TNumNodes, class TWallModel >
400 inline std::ostream& operator << (std::ostream& rOStream, const TwoFluidNavierStokesWallCondition<TDim,TNumNodes,TWallModel>& rThis)
401 {
402  rThis.PrintInfo(rOStream);
403  rOStream << std::endl;
404  rThis.PrintData(rOStream);
405 
406  return rOStream;
407 }
408 
411 } // namespace Kratos.
Condition & operator=(Condition const &rOther)
Assignment operator.
Definition: condition.h:181
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
Geometry base class.
Definition: geometry.h:71
Implements a wall condition for the Navier-Stokes (and Stokes) monolithic formulations This condition...
Definition: navier_stokes_wall_condition.h:72
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
Properties encapsulates data shared by different Elements or Conditions. It can store any type of dat...
Definition: properties.h:69
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
Implements a wall condition for the Navier-Stokes monolithic formulation.
Definition: two_fluid_navier_stokes_wall_condition.h:66
Condition::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
Create a new TwoFluidNavierStokesWallCondition object.
Definition: two_fluid_navier_stokes_wall_condition.h:182
Vector VectorType
Definition: two_fluid_navier_stokes_wall_condition.h:88
TwoFluidNavierStokesWallCondition(IndexType NewId=0)
Default constructor.
Definition: two_fluid_navier_stokes_wall_condition.h:108
Condition::Pointer Create(IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties) const override
Create a new TwoFluidNavierStokesWallCondition object.
Definition: two_fluid_navier_stokes_wall_condition.h:196
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(TwoFluidNavierStokesWallCondition)
Pointer definition of TwoFluidNavierStokesWallCondition.
std::vector< std::size_t > EquationIdVectorType
Definition: two_fluid_navier_stokes_wall_condition.h:94
std::string Info() const override
Turn back information as a string.
Definition: two_fluid_navier_stokes_wall_condition.h:245
TwoFluidNavierStokesWallCondition(TwoFluidNavierStokesWallCondition const &rOther)
Copy constructor.
Definition: two_fluid_navier_stokes_wall_condition.h:152
Geometry< NodeType >::PointsArrayType NodesArrayType
Definition: two_fluid_navier_stokes_wall_condition.h:86
NavierStokesWallCondition< TDim, TNumNodes, TWallModel... > BaseType
Definition: two_fluid_navier_stokes_wall_condition.h:74
TwoFluidNavierStokesWallCondition(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
Constructor using Properties.
Definition: two_fluid_navier_stokes_wall_condition.h:143
std::vector< Dof< double >::Pointer > DofsVectorType
Definition: two_fluid_navier_stokes_wall_condition.h:96
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: two_fluid_navier_stokes_wall_condition.h:259
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: two_fluid_navier_stokes_wall_condition.h:253
std::size_t IndexType
Definition: two_fluid_navier_stokes_wall_condition.h:92
~TwoFluidNavierStokesWallCondition() override
Destructor.
Definition: two_fluid_navier_stokes_wall_condition.h:158
Condition::Pointer Clone(IndexType NewId, NodesArrayType const &rThisNodes) const override
Definition: two_fluid_navier_stokes_wall_condition.h:210
TwoFluidNavierStokesWallCondition & operator=(TwoFluidNavierStokesWallCondition const &rOther)
Assignment operator.
Definition: two_fluid_navier_stokes_wall_condition.h:166
Node NodeType
Definition: two_fluid_navier_stokes_wall_condition.h:80
Geometry< NodeType > GeometryType
Definition: two_fluid_navier_stokes_wall_condition.h:84
TwoFluidNavierStokesWallCondition(IndexType NewId, GeometryType::Pointer pGeometry)
Constructor using Geometry.
Definition: two_fluid_navier_stokes_wall_condition.h:130
Matrix MatrixType
Definition: two_fluid_navier_stokes_wall_condition.h:90
BaseType::ConditionDataStruct ConditionDataStruct
Definition: two_fluid_navier_stokes_wall_condition.h:78
TwoFluidNavierStokesWallCondition(IndexType NewId, const NodesArrayType &ThisNodes)
Constructor using an array of nodes.
Definition: two_fluid_navier_stokes_wall_condition.h:118
Properties PropertiesType
Definition: two_fluid_navier_stokes_wall_condition.h:82
#define KRATOS_SERIALIZE_SAVE_BASE_CLASS(Serializer, BaseType)
Definition: define.h:812
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
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
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