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.
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 #include "includes/define.h"
25 #include "includes/condition.h"
26 #include "includes/model_part.h"
27 #include "includes/serializer.h"
28 #include "includes/process_info.h"
29 
30 // Application includes
32 
33 namespace Kratos
34 {
35 
38 
41 
45 
49 
53 
57 
70 template<unsigned int TDim, unsigned int TNumNodes, class... TWallModel>
71 class KRATOS_API(FLUID_DYNAMICS_APPLICATION) NavierStokesWallCondition : public Condition
72 {
73 public:
76 
79 
81  {
82  double wGauss; // Gauss point weight
83  array_1d<double, 3> Normal; // Condition normal
84  array_1d<double, TNumNodes> N; // Gauss point shape functions values
85  Vector ViscousStress; // Viscous stresses that are retrieved from parent
86  };
87 
88  static constexpr std::size_t VoigtSize = 3 * (TDim-1);
89  static constexpr std::size_t BlockSize = TDim + 1;
90  static constexpr std::size_t LocalSize = TNumNodes*BlockSize;
91 
92  using Condition::SizeType;
93 
94  typedef Node NodeType;
95 
97 
99 
101 
103 
105 
106  typedef std::size_t IndexType;
107 
108  typedef std::vector<std::size_t> EquationIdVectorType;
109 
110  typedef std::vector< Dof<double>::Pointer > DofsVectorType;
111 
115 
117 
123  {
124  }
125 
127 
132  Condition(NewId,ThisNodes)
133  {
134  }
135 
137 
141  NavierStokesWallCondition(IndexType NewId, GeometryType::Pointer pGeometry):
142  Condition(NewId,pGeometry)
143  {
144  }
145 
147 
152  NavierStokesWallCondition(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties):
153  Condition(NewId,pGeometry,pProperties)
154  {
155  }
156 
159  Condition(rOther)
160  {
161  }
162 
165 
166 
170 
173  {
174  Condition::operator=(rOther);
175  return *this;
176  }
177 
181 
183 
188  Condition::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes, PropertiesType::Pointer pProperties) const override
189  {
190  return Kratos::make_intrusive<NavierStokesWallCondition>(NewId, GetGeometry().Create(ThisNodes), pProperties);
191  }
192 
194 
199  Condition::Pointer Create(IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties) const override
200  {
201  return Kratos::make_intrusive< NavierStokesWallCondition >(NewId, pGeom, pProperties);
202  }
203 
210  Condition::Pointer Clone(IndexType NewId, NodesArrayType const& rThisNodes) const override
211  {
212  Condition::Pointer pNewCondition = Create(NewId, GetGeometry().Create( rThisNodes ), pGetProperties() );
213 
214  pNewCondition->SetData(this->GetData());
215  pNewCondition->SetFlags(this->GetFlags());
216 
217  return pNewCondition;
218  }
219 
220 
222 
228  void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix,
229  VectorType& rRightHandSideVector,
230  const ProcessInfo& rCurrentProcessInfo) override;
231 
232 
234 
239  void CalculateLeftHandSide(MatrixType& rLeftHandSideMatrix,
240  const ProcessInfo& rCurrentProcessInfo) override;
241 
242 
244 
249  void CalculateRightHandSide(VectorType& rRightHandSideVector,
250  const ProcessInfo& rCurrentProcessInfo) override;
251 
252 
253 
255 
258  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
259 
266  void EquationIdVector(
267  EquationIdVectorType& rResult,
268  const ProcessInfo& rCurrentProcessInfo) const override;
269 
276  void GetDofList(
277  DofsVectorType& rConditionDofList,
278  const ProcessInfo& rCurrentProcessInfo) const override;
279 
280  void Calculate(
281  const Variable< array_1d<double,3> >& rVariable,
282  array_1d<double,3>& Output,
283  const ProcessInfo& rCurrentProcessInfo) override;
284 
288 
289 
293 
294 
298 
300  std::string Info() const override
301  {
302  std::stringstream buffer;
303  buffer << "NavierStokesWallCondition" << TDim << "D";
304  return buffer.str();
305  }
306 
308  void PrintInfo(std::ostream& rOStream) const override
309  {
310  rOStream << "NavierStokesWallCondition";
311  }
312 
314  void PrintData(std::ostream& rOStream) const override {}
315 
316 
320 
321 
323 
324 protected:
327 
328 
332 
333 
337 
338 
342 
348  void CalculateNormal(array_1d<double,3>& rAreaNormal);
349 
359  void ComputeGaussPointLHSContribution(
361  const ConditionDataStruct& rData,
362  const ProcessInfo& rProcessInfo);
363 
373  void ComputeGaussPointRHSContribution(
375  const ConditionDataStruct& rData,
376  const ProcessInfo& rProcessInfo);
377 
386  void ComputeRHSNeumannContribution(
388  const ConditionDataStruct& data);
389 
399  void ComputeRHSOutletInflowContribution(
401  const ConditionDataStruct& rData,
402  const ProcessInfo& rProcessInfo);
403 
407 
408 
412 
413 
417 
418 
420 
421 private:
424 
425 
429 
430 
434 
435  friend class Serializer;
436 
437  void save(Serializer& rSerializer) const override
438  {
440  }
441 
442  void load(Serializer& rSerializer) override
443  {
445  }
446 
450 
451 
455 
464  void CalculateGaussPointSlipTangentialCorrectionLHSContribution(
465  BoundedMatrix<double,LocalSize,LocalSize>& rLeftHandSideMatrix,
466  const ConditionDataStruct& rDataStruct);
467 
476  void CalculateGaussPointSlipTangentialCorrectionRHSContribution(
477  array_1d<double,LocalSize>& rRightHandSideVector,
478  const ConditionDataStruct& rDataStruct);
479 
488  void ProjectViscousStress(
489  const Vector& rViscousStress,
490  const array_1d<double,3> rNormal,
491  array_1d<double,3>& rProjectedViscousStress);
492 
499  void SetTangentialProjectionMatrix(
500  const array_1d<double,3>& rUnitNormal,
501  BoundedMatrix<double,TDim,TDim>& rTangProjMat)
502  {
503  noalias(rTangProjMat) = IdentityMatrix(TDim,TDim);
504  for (std::size_t d1 = 0; d1 < TDim; ++d1) {
505  for (std::size_t d2 = 0; d2 < TDim; ++d2) {
506  rTangProjMat(d1,d2) -= rUnitNormal[d1]*rUnitNormal[d2];
507  }
508  }
509  }
510 
511  template<typename TWallModelType>
512  int WallModelCheckCall(const ProcessInfo& rProcessInfo) const
513  {
514  return TWallModelType::Check(this, rProcessInfo);
515  }
516 
517  template<typename TWallModelType>
518  void AddWallModelRightHandSideCall(
519  VectorType& rRHS,
520  const ProcessInfo& rProcessInfo)
521  {
522  TWallModelType::AddWallModelRightHandSide(rRHS, this, rProcessInfo);
523  }
524 
525  template<typename TWallModelType>
526  void AddWallModelLeftHandSideCall(
527  MatrixType& rLHS,
528  const ProcessInfo& rProcessInfo)
529  {
530  TWallModelType::AddWallModelLeftHandSide(rLHS, this, rProcessInfo);
531  }
532 
533  template<typename TWallModelType>
534  void AddWallModelLocalSystemCall(
535  MatrixType& rLHS,
536  VectorType& rRHS,
537  const ProcessInfo& rProcessInfo)
538  {
539  TWallModelType::AddWallModelLocalSystem(rLHS, rRHS, this, rProcessInfo);
540  }
541 
545 
546 
550 
551 
555 
556 
558 
559 }; // Class NavierStokesWallCondition
560 
561 
563 
566 
567 
571 
572 
574 template< unsigned int TDim, unsigned int TNumNodes, class TWallModel >
575 inline std::istream& operator >> (std::istream& rIStream, NavierStokesWallCondition<TDim,TNumNodes,TWallModel>& rThis)
576 {
577  return rIStream;
578 }
579 
581 template< unsigned int TDim, unsigned int TNumNodes, class TWallModel >
582 inline std::ostream& operator << (std::ostream& rOStream, const NavierStokesWallCondition<TDim,TNumNodes,TWallModel>& rThis)
583 {
584  rThis.PrintInfo(rOStream);
585  rOStream << std::endl;
586  rThis.PrintData(rOStream);
587 
588  return rOStream;
589 }
590 
592 
594 
595 
596 } // namespace Kratos.
Base class for all Conditions.
Definition: condition.h:59
std::size_t SizeType
Definition: condition.h:94
Condition & operator=(Condition const &rOther)
Assignment operator.
Definition: condition.h:181
std::size_t IndexType
Definition: flags.h:74
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
Vector VectorType
Definition: navier_stokes_wall_condition.h:102
~NavierStokesWallCondition() override
Destructor.
Definition: navier_stokes_wall_condition.h:164
Condition::Pointer Clone(IndexType NewId, NodesArrayType const &rThisNodes) const override
Definition: navier_stokes_wall_condition.h:210
NavierStokesWallCondition(IndexType NewId=0)
Default constructor.
Definition: navier_stokes_wall_condition.h:122
std::size_t IndexType
Definition: navier_stokes_wall_condition.h:106
Matrix MatrixType
Definition: navier_stokes_wall_condition.h:104
std::string Info() const override
Turn back information as a string.
Definition: navier_stokes_wall_condition.h:300
Condition::Pointer Create(IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties) const override
Create a new NavierStokesWallCondition object.
Definition: navier_stokes_wall_condition.h:199
NavierStokesWallCondition & operator=(NavierStokesWallCondition const &rOther)
Assignment operator.
Definition: navier_stokes_wall_condition.h:172
NavierStokesWallCondition(IndexType NewId, GeometryType::Pointer pGeometry)
Constructor using Geometry.
Definition: navier_stokes_wall_condition.h:141
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(NavierStokesWallCondition)
Pointer definition of NavierStokesWallCondition.
std::vector< std::size_t > EquationIdVectorType
Definition: navier_stokes_wall_condition.h:108
Node NodeType
Definition: navier_stokes_wall_condition.h:94
Properties PropertiesType
Definition: navier_stokes_wall_condition.h:96
NavierStokesWallCondition(NavierStokesWallCondition const &rOther)
Copy constructor.
Definition: navier_stokes_wall_condition.h:158
Geometry< NodeType >::PointsArrayType NodesArrayType
Definition: navier_stokes_wall_condition.h:100
std::vector< Dof< double >::Pointer > DofsVectorType
Definition: navier_stokes_wall_condition.h:110
NavierStokesWallCondition(IndexType NewId, const NodesArrayType &ThisNodes)
Constructor using an array of nodes.
Definition: navier_stokes_wall_condition.h:131
Geometry< NodeType > GeometryType
Definition: navier_stokes_wall_condition.h:98
NavierStokesWallCondition(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
Constructor using Properties.
Definition: navier_stokes_wall_condition.h:152
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: navier_stokes_wall_condition.h:314
Condition::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
Create a new NavierStokesWallCondition object.
Definition: navier_stokes_wall_condition.h:188
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: navier_stokes_wall_condition.h:308
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
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
#define KRATOS_SERIALIZE_SAVE_BASE_CLASS(Serializer, BaseType)
Definition: define.h:812
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
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
TDataType Calculate(GeometryType &dummy, const Variable< TDataType > &rVariable)
Definition: add_geometries_to_python.cpp:103
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
AMatrix::IdentityMatrix< double > IdentityMatrix
Definition: amatrix_interface.h:564
Internals::Matrix< double, AMatrix::dynamic, 1 > Vector
Definition: amatrix_interface.h:472
array_1d< double, 3 > VectorType
Definition: solid_mechanics_application_variables.cpp:19
std::istream & operator>>(std::istream &rIStream, LinearMasterSlaveConstraint &rThis)
input stream function
T & noalias(T &TheMatrix)
Definition: amatrix_interface.h:484
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
ProcessInfo
Definition: edgebased_PureConvection.py:116
data
Definition: mesh_to_mdpa_converter.py:59
def load(f)
Definition: ode_solve.py:307
subroutine d1(DSTRAN, D, dtime, NDI, NSHR, NTENS)
Definition: TensorModule.f:594
Definition: navier_stokes_wall_condition.h:81
array_1d< double, TNumNodes > N
Definition: navier_stokes_wall_condition.h:84
Vector ViscousStress
Definition: navier_stokes_wall_condition.h:85
double wGauss
Definition: navier_stokes_wall_condition.h:82
array_1d< double, 3 > Normal
Definition: navier_stokes_wall_condition.h:83