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.
stationary_stokes.h
Go to the documentation of this file.
1 #ifndef KRATOS_STATIONARY_STOKES_H
2 #define KRATOS_STATIONARY_STOKES_H
3 
4 /*
5 ==============================================================================
6 Kratos
7 A General Purpose Software for Multi-Physics Finite Element Analysis
8 Version 1.0 (Released on march 05, 2007).
9 
10 Copyright 2007
11 Pooyan Dadvand, Riccardo Rossi
12 pooyan@cimne.upc.edu
13 rrossi@cimne.upc.edu
14 CIMNE (International Center for Numerical Methods in Engineering),
15 Gran Capita' s/n, 08034 Barcelona, Spain
16 
17 Permission is hereby granted, free of charge, to any person obtaining
18 a copy of this software and associated documentation files (the
19 "Software"), to deal in the Software without restriction, including
20 without limitation the rights to use, copy, modify, merge, publish,
21 distribute, sublicense and/or sell copies of the Software, and to
22 permit persons to whom the Software is furnished to do so, subject to
23 the following condition:
24 
25 Distribution of this code for any commercial purpose is permissible
26 ONLY BY DIRECT ARRANGEMENT WITH THE COPYRIGHT OWNER.
27 
28 The above copyright notice and this permission notice shall be
29 included in all copies or substantial portions of the Software.
30 
31 THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND,
32 EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF
33 MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT.
34 IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY
35 CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT,
36 TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE
37 SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.
38 
39 ==============================================================================
40  */
41 
42 #ifndef KRATOS_FLUID_TAYLOR_HOOD_H
43 #define KRATOS_FLUID_TAYLOR_HOOD_H
44 // System includes
45 #include <string>
46 #include <iostream>
47 
48 
49 // External includes
50 
51 
52 // Project includes
53 #include "includes/define.h"
54 #include "includes/serializer.h"
55 #include "includes/element.h"
59 
60 
61 namespace Kratos
62 {
65 
68 
72 
76 
80 
84 
86 template< unsigned int TDim >
87 class StationaryStokes : public Element
88 {
89 public:
92 
95 
98 
101 
104 
108 
111  Element(NewId)
112  {}
113 
115  StationaryStokes(IndexType NewId, GeometryType::Pointer pGeometry) :
116  Element(NewId, pGeometry)
117  {}
118 
120  StationaryStokes(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties) :
121  Element(NewId, pGeometry, pProperties)
122  {}
123 
125  ~StationaryStokes() override
126  {}
127 
131 
132 
136 
138  Element::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes, PropertiesType::Pointer pProperties) const override;
139 
140 
142 
147  Element::Pointer Create(IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties) const override;
148 
150 
152  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
153 
155  void Initialize(const ProcessInfo &rCurrentProcessInfo) override;
156 
158 
163  void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix, VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo) override;
164 
165  void CalculateRightHandSide(VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo) override
166  {
167  MatrixType TmpLHS;
168  this->CalculateLocalSystem(TmpLHS,rRightHandSideVector,rCurrentProcessInfo);
169  }
170 
172  void GetDofList(DofsVectorType& rElementalDofList, const ProcessInfo& rCurrentProcessInfo) const override;
173 
175  void EquationIdVector(Element::EquationIdVectorType& rResult, const ProcessInfo& rCurrentProcessInfo) const override;
176 
177  void GetFirstDerivativesVector(Vector &rValues, int Step = 0) const override;
178 
182 
183 
187 
188 
192 
194  std::string Info() const override
195  {
196  std::stringstream buffer;
197  buffer << "StationaryStokes" << this->GetGeometry().WorkingSpaceDimension() << "D #" << Id();
198  return buffer.str();
199  }
200 
202  void PrintInfo(std::ostream& rOStream) const override
203  {
204  rOStream << "StationaryStokes" << this->GetGeometry().WorkingSpaceDimension() << "D #" << Id() << std::endl;
205  rOStream << "Number of Nodes: " << this->GetGeometry().PointsNumber() << std::endl;
206  rOStream << "Integration method: " << static_cast<int>(this->mIntegrationMethod);
207  }
208 
210  void PrintData(std::ostream& rOStream) const override
211  {
212  this->PrintInfo(rOStream);
213  rOStream << "Geometry Data: " << std::endl;
214  this->GetGeometry().PrintData(rOStream);
215  }
216 
217 
221 
222 
224 
225 protected:
228 
229 
233 
234 
238 
239 
243 
244  void AddMomentumTerms(MatrixType &rLHS,
245  VectorType &rRHS,
246  const double Density,
247  const double Viscosity,
248  const array_1d<double,3>& BodyForce,
249  const double TauTwo,
250  const ShapeFunctionsType &N,
251  const ShapeDerivativesType &DN_DX,
252  const double Weigth);
253 
254  void AddContinuityTerms(MatrixType &rLHS,
255  VectorType &rRHS,
256  const double Density,
257  const array_1d<double,3>& BodyForce,
258  const double TauOne,
259  const ShapeFunctionsType &N,
260  const ShapeDerivativesType &DN_DX,
261  const double Weight);
262 
263  template< class TVariableType >
264  void EvaluateInPoint(TVariableType& rResult,
266  const ShapeFunctionsType& rShapeFunc,
267  GeometryType& rGeom)
268  {
269  rResult = rShapeFunc[0] * rGeom[0].FastGetSolutionStepValue(Var);
270 
271  for(SizeType i = 1; i < rShapeFunc.size(); i++)
272  {
273  rResult += rShapeFunc[i] * rGeom[i].FastGetSolutionStepValue(Var);
274  }
275  }
276 
277  virtual void CalculateTau(double& TauOne,
278  double& TauTwo,
279  const double DynViscosity);
280 
281  double ElementSize();
282 
286 
287 
291 
292 
296 
297 
299 
300 private:
303 
304 
308 
310  Element::IntegrationMethod mIntegrationMethod;
311 
313  std::vector< ShapeDerivativesType > mDN_DX;
314 
316  std::vector< double > mGaussWeight;
317 
321 
322  friend class Serializer;
323 
324  void save(Serializer& rSerializer) const override
325  {
326  KRATOS_TRY;
328 
329  unsigned int IntMethod = 0;
330  switch(mIntegrationMethod)
331  {
333  IntMethod = 1;
334  break;
336  IntMethod = 2;
337  break;
339  IntMethod = 3;
340  break;
342  IntMethod = 4;
343  break;
345  IntMethod = 5;
346  break;
347  default:
348  KRATOS_ERROR << "Unknown integration method encountered on serializer save for StationaryStokes element: " << static_cast<int>(mIntegrationMethod) << std::endl;
349  break;
350  }
351  rSerializer.save("IntMethod",IntMethod);
352  rSerializer.save("mDN_DX",mDN_DX);
353  rSerializer.save("mGaussWeight",mGaussWeight);
354  KRATOS_CATCH("");
355  }
356 
357  void load(Serializer& rSerializer) override
358  {
359  KRATOS_TRY;
361 
362  unsigned int IntMethod = 0;
363  rSerializer.load("IntMethod",IntMethod);
364  switch(static_cast<int>(mIntegrationMethod))
365  {
366  case 1:
367  mIntegrationMethod = GeometryData::IntegrationMethod::GI_GAUSS_1;
368  break;
369  case 2:
370  mIntegrationMethod = GeometryData::IntegrationMethod::GI_GAUSS_2;
371  break;
372  case 3:
373  mIntegrationMethod = GeometryData::IntegrationMethod::GI_GAUSS_3;
374  break;
375  case 4:
376  mIntegrationMethod = GeometryData::IntegrationMethod::GI_GAUSS_4;
377  break;
378  case 5:
379  mIntegrationMethod = GeometryData::IntegrationMethod::GI_GAUSS_5;
380  break;
381  default:
382  KRATOS_ERROR << "Unknown integration method encountered on serializer load for StationaryStokes element: " << static_cast<int>(IntMethod);
383  break;
384  }
385  rSerializer.load("mDN_DX",mDN_DX);
386  rSerializer.load("mGaussWeight",mGaussWeight);
387  KRATOS_CATCH("");
388  }
389 
393 
394 
398 
399 
403 
404 
408 
409 
413 
416 
418  StationaryStokes(StationaryStokes const& rOther);
419 
420 
422 
423 }; // Class StationaryStokes
424 
426 
429 
430 
434 
435 
437 template< unsigned int TDim >
438 inline std::istream& operator >> (std::istream& rIStream,
439  StationaryStokes<TDim>& rThis)
440 {
441  return rIStream;
442 }
443 
445 template< unsigned int TDim >
446 inline std::ostream& operator << (std::ostream& rOStream,
447  const StationaryStokes<TDim>& rThis)
448 {
449  rThis.PrintInfo(rOStream);
450  rOStream << std::endl;
451  rThis.PrintData(rOStream);
452 
453  return rOStream;
454 }
456 
458 
459 } // namespace Kratos.
460 
461 #endif // KRATOS_FLUID_TAYLOR_HOOD_H
462 
463 
464 #endif // KRATOS_STATIONARY_STOKES_H
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
Base class for all Elements.
Definition: element.h:60
std::size_t SizeType
Definition: element.h:94
Element(IndexType NewId=0)
Definition: element.h:121
std::vector< DofType::Pointer > DofsVectorType
Definition: element.h:100
Matrix MatrixType
Definition: element.h:90
std::vector< std::size_t > EquationIdVectorType
Definition: element.h:98
std::size_t IndexType
Definition: flags.h:74
GeometryType & GetGeometry()
Returns the reference of the geometry.
Definition: geometrical_object.h:158
IntegrationMethod
Definition: geometry_data.h:76
Geometry base class.
Definition: geometry.h:71
SizeType PointsNumber() const
Definition: geometry.h:528
virtual void PrintData(std::ostream &rOStream) const
Print object's data.
Definition: geometry.h:3834
SizeType WorkingSpaceDimension() const
Definition: geometry.h:1287
IndexType Id() const
Definition: indexed_object.h:107
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
void save(std::string const &rTag, std::array< TDataType, TDataSize > const &rObject)
Definition: serializer.h:545
GLS-stabilized element for the solution of the stationary Stokes problem.
Definition: stationary_stokes.h:88
void EquationIdVector(Element::EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const override
Fill given vector with the linear system row index for the element's degrees of freedom.
Definition: stationary_stokes.cpp:179
Element::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
Create a new StationaryStokes element and return a pointer to it.
Definition: stationary_stokes.cpp:8
int Check(const ProcessInfo &rCurrentProcessInfo) const override
Check that all required data containers are properly initialized and registered in Kratos.
Definition: stationary_stokes.cpp:20
void GetFirstDerivativesVector(Vector &rValues, int Step=0) const override
Definition: stationary_stokes.cpp:201
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: stationary_stokes.h:210
Kratos::Geometry< Node > GeometryType
Type for geometry calls.
Definition: stationary_stokes.h:103
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(StationaryStokes)
Pointer definition of StationaryStokes.
StationaryStokes(IndexType NewId, GeometryType::Pointer pGeometry)
Constructor using a Geometry instance.
Definition: stationary_stokes.h:115
void EvaluateInPoint(TVariableType &rResult, const Kratos::Variable< TVariableType > &Var, const ShapeFunctionsType &rShapeFunc, GeometryType &rGeom)
Definition: stationary_stokes.h:264
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: stationary_stokes.h:202
StationaryStokes(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
Constructor using geometry and properties.
Definition: stationary_stokes.h:120
void AddContinuityTerms(MatrixType &rLHS, VectorType &rRHS, const double Density, const array_1d< double, 3 > &BodyForce, const double TauOne, const ShapeFunctionsType &N, const ShapeDerivativesType &DN_DX, const double Weight)
Definition: stationary_stokes.cpp:274
StationaryStokes(IndexType NewId=0)
Default constructor.
Definition: stationary_stokes.h:110
Kratos::Vector ShapeFunctionsType
Type for shape function values container.
Definition: stationary_stokes.h:97
void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Evaluate the elemental contribution to the problem.
Definition: stationary_stokes.cpp:96
Kratos::Matrix ShapeDerivativesType
Type for shape function derivatives container.
Definition: stationary_stokes.h:100
void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: stationary_stokes.h:165
void AddMomentumTerms(MatrixType &rLHS, VectorType &rRHS, const double Density, const double Viscosity, const array_1d< double, 3 > &BodyForce, const double TauTwo, const ShapeFunctionsType &N, const ShapeDerivativesType &DN_DX, const double Weigth)
Definition: stationary_stokes.cpp:222
Element::Pointer Create(IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties) const override
Create a new element of this type.
void GetDofList(DofsVectorType &rElementalDofList, const ProcessInfo &rCurrentProcessInfo) const override
Fill given array with containing the element's degrees of freedom.
Definition: stationary_stokes.cpp:156
friend class Serializer
Definition: stationary_stokes.h:322
~StationaryStokes() override
Destructor.
Definition: stationary_stokes.h:125
void Initialize(const ProcessInfo &rCurrentProcessInfo) override
Calculate Shape function derivatives and Jacobian at each integration point.
Definition: stationary_stokes.cpp:54
virtual void CalculateTau(double &TauOne, double &TauTwo, const double DynViscosity)
Definition: stationary_stokes.cpp:321
std::string Info() const override
Turn back information as a string.
Definition: stationary_stokes.h:194
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_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
std::size_t IndexType
The definition of the index type.
Definition: key_hash.h:35
#define KRATOS_ERROR
Definition: exception.h:161
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
ModelPart::NodesContainerType NodesArrayType
Definition: gid_gauss_point_container.h:42
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
N
Definition: sensitivityMatrix.py:29
integer i
Definition: TensorModule.f:17