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.
dynamic_vms.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: Jordi Cotela
11 //
12 
13 
14 #if !defined(KRATOS_DYNAMIC_VMS_H_INCLUDED )
15 #define KRATOS_DYNAMIC_VMS_H_INCLUDED
16 
17 // System includes
18 #include <string>
19 #include <iostream>
20 
21 
22 // External includes
23 
24 
25 // Project includes
26 #include "containers/array_1d.h"
27 #include "includes/define.h"
28 #include "includes/element.h"
29 #include "includes/serializer.h"
30 #include "geometries/geometry.h"
31 #include "utilities/math_utils.h"
32 
33 // Application includes
35 
36 namespace Kratos
37 {
38 
41 
44 
48 
52 
56 
60 
62 
104 template< unsigned int TDim >
105 class DynamicVMS : public Element
106 {
107 public:
110 
113 
116 
119 
123 
124  //Constructors.
125 
127 
131  DynamicVMS(IndexType NewId, const NodesArrayType& ThisNodes);
132 
134 
138  DynamicVMS(IndexType NewId, GeometryType::Pointer pGeometry);
139 
141 
147  DynamicVMS(IndexType NewId, GeometryType::Pointer pGeometry, const GeometryData::IntegrationMethod ThisIntegrationMethod);
148 
150 
155  DynamicVMS(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties);
156 
165  DynamicVMS(IndexType NewId, GeometryType::Pointer pGeometry, Properties::Pointer pProperties, const GeometryData::IntegrationMethod ThisIntegrationMethod);
166 
168  ~DynamicVMS() override;
169 
173 
174 
178 
180 
188  Element::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes,
189  PropertiesType::Pointer pProperties) const override;
190 
192 
197  Element::Pointer Create(
198  IndexType NewId,
199  GeometryType::Pointer pGeom,
200  PropertiesType::Pointer pProperties) const override;
201 
203  void Initialize(const ProcessInfo &rCurrentProcessInfo) override;
204 
210  void InitializeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
211 
213 
216  void InitializeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo) override;
217 
218 
219  void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix,
220  VectorType &rRightHandSideVector,
221  const ProcessInfo &rCurrentProcessInfo) override;
222 
223 
224  void CalculateRightHandSide(VectorType& rRightHandSideVector,
225  const ProcessInfo& rCurrentProcessInfo) override;
226 
227 
229 
237  void CalculateLocalVelocityContribution(MatrixType& rDampingMatrix,
238  VectorType& rRightHandSideVector,
239  const ProcessInfo& rCurrentProcessInfo) override;
240 
241  void CalculateMassMatrix(MatrixType &rMassMatrix, const ProcessInfo &rCurrentProcessInfo) override;
242 
243 
245 
252  const ProcessInfo& rCurrentProcessInfo) const override;
253 
255 
259  void GetDofList(DofsVectorType& rElementalDofList,
260  const ProcessInfo& rCurrentProcessInfo) const override;
261 
263 
267  void GetFirstDerivativesVector(Vector& rValues, int Step = 0) const override;
268 
270 
274  void GetSecondDerivativesVector(Vector& rValues, int Step = 0) const override;
275 
276 
277 // virtual void Calculate(const Variable<double>& rVariable,
278 // double& rOutput,
279 // const ProcessInfo& rCurrentProcessInfo);
280 
281  void Calculate(const Variable< array_1d<double,3> >& rVariable,
282  array_1d<double,3>& rOutput,
283  const ProcessInfo& rCurrentProcessInfo) override;
284 
285 
287  std::vector<array_1d<double, 3 > >& rOutput,
288  const ProcessInfo& rCurrentProcessInfo) override;
289 
290  void CalculateOnIntegrationPoints(const Variable<double>& rVariable,
291  std::vector<double>& rValues,
292  const ProcessInfo& rCurrentProcessInfo) override;
293 
296  std::vector<array_1d<double, 6 > >& rValues,
297  const ProcessInfo& rCurrentProcessInfo) override
298  {}
299 
302  std::vector<Vector>& rValues,
303  const ProcessInfo& rCurrentProcessInfo) override
304  {}
305 
308  std::vector<Matrix>& rValues,
309  const ProcessInfo& rCurrentProcessInfo) override
310  {}
311 
312  void SetValuesOnIntegrationPoints(const Variable<double> &rVariable, const std::vector<double> &rValues, const ProcessInfo &rCurrentProcessInfo) override;
313 
317 
319 
322 
326 
327 
331 
333  std::string Info() const override;
334 
336  void PrintInfo(std::ostream& rOStream) const override;
337 
338  void PrintData(std::ostream &rOStream) const override;
339 
343 
344 
346 
347 protected:
350 
351 
355 
356 
360 
361 
365 
367  virtual void CalculateGeometryData();
368 
370 
373  virtual void UpdateSubscale(const ProcessInfo& rCurrentProcessInfo);
374 
375  virtual void LinearUpdateSubscale(const ProcessInfo& rCurrentProcessInfo);
376 
377  virtual void CalculateASGSVelocityContribution(MatrixType& rDampingMatrix,
378  VectorType& rRightHandSideVector,
379  const ProcessInfo& rCurrentProcessInfo);
380 
381  virtual void CalculateOSSVelocityContribution(MatrixType& rDampingMatrix,
382  VectorType& rRightHandSideVector,
383  const ProcessInfo& rCurrentProcessInfo);
384 
385  virtual void LumpedMassMatrix(MatrixType& rMassMatrix);
386 
387  virtual void ConsistentMassMatrix(MatrixType& rMassMatrix);
388 
389 
390  template< typename TValueType >
391  void EvaluateInPoint(TValueType& rValue,
392  const Kratos::Variable< TValueType >& rVariable,
393  const ShapeFunctionsType& rN)
394  {
395  GeometryType& rGeom = this->GetGeometry();
396  const unsigned int NumNodes = rGeom.PointsNumber();
397 
398  rValue = rN[0] * rGeom[0].FastGetSolutionStepValue(rVariable);
399  for (unsigned int i = 1; i < NumNodes; i++)
400  rValue += rN[i] * rGeom[i].FastGetSolutionStepValue(rVariable);
401  }
402 
403  template< typename TValueType >
404  void EvaluateInPoint(TValueType& rValue,
405  const Kratos::Variable< TValueType >& rVariable,
406  const ShapeFunctionsType& rN,
407  const unsigned int Step)
408  {
409  const GeometryType& rGeom = this->GetGeometry();
410  const unsigned int NumNodes = rGeom.PointsNumber();
411 
412  rValue = rN[0] * rGeom[0].FastGetSolutionStepValue(rVariable,Step);
413  for (unsigned int i = 1; i < NumNodes; i++)
414  rValue += rN[i] * rGeom[i].FastGetSolutionStepValue(rVariable,Step);
415  }
416 
417  virtual void EvaluateConvVelocity(array_1d<double,3>& rConvection,
418  const ShapeFunctionsType& rN);
419 
420  virtual void EvaluateConvVelocity(array_1d<double,3>& rConvection,
421  const array_1d<double,3>& rSubscaleVel,
422  const ShapeFunctionsType& rN);
423 
430  virtual void EvaluateViscosity(double& rViscosity,
431  const ShapeFunctionsType& rN);
432 
433  virtual void ConvectionOperator(Vector& rResult,
434  const array_1d<double,3>& rConvVel);
435 
437 
443  virtual double TauOne(const double Density,
444  const double Viscosity,
445  const double ConvVel);
446 
448 
454  virtual double TauTwo(const double Density,
455  const double Viscosity,
456  const double ConvVel);
457 
459 
466  virtual double TauTime(const double Density,
467  const double Viscosity,
468  const double ConvVel,
469  const double Dt);
470 
472  virtual double InvTauTime(const double Density,
473  const double Viscosity,
474  const double ConvVel,
475  const double Dt);
476 
477 
478  virtual void ASGSMomentumResidual(array_1d<double,3>& rResult,
479  const double Density,
480  const array_1d<double,3>& rConvVel,
481  const ShapeFunctionsType& rN);
482 
483  virtual void OSSMomentumResidual(array_1d<double,3>& rResult,
484  const double Density,
485  const array_1d<double,3>& rConvVel,
486  const ShapeFunctionsType& rN);
487 
488  virtual void MassResidual(double& rResult);
489 
490  virtual void AddViscousTerm(MatrixType& rDampingMatrix,
491  const double Weight,
492  const ShapeDerivativesType& rDN_DX);
493 
494  virtual void DenseSystemSolve(const Matrix& rA,
495  const Vector& rB,
496  Vector& rX);
497 
499  const ShapeDerivativesType& rDN_DX);
500 
501 
505 
506 
510 
511 
515 
516 
518 
519 private:
522 
524  static const double mSubscaleTol;
525 
527  static const double mSubscaleRHSTol;
528 
532 
533  GeometryData::IntegrationMethod mIntegrationMethod;
534 
535  ShapeDerivativesType mDN_DX;
536 
537  double mDetJ;
538 
539  double mElemSize;
540 
542  std::vector< array_1d<double,3> > mSubscaleVel;
543 
545  std::vector< array_1d<double,3> > mOldSubscaleVel;
546 
548  std::vector< unsigned int > mIterCount;
549 
550 // std::vector< array_1d<double,3> > mExtraTerm;
551 
555 
556  friend class Serializer;
557 
559  DynamicVMS();
560 
561  void save(Serializer& rSerializer) const override;
562 
563 
564  void load(Serializer& rSerializer) override;
565 
566 
570 
571 
575 
576 
580 
581 
585 
586 
590 
592  DynamicVMS & operator=(DynamicVMS const& rOther);
593 
595  DynamicVMS(DynamicVMS const& rOther);
596 
598 
599 }; // Class DynamicVMS
600 
602 
605 
606 
610 
611 
613 template< unsigned int TDim >
614 inline std::istream& operator >>(std::istream& rIStream,
615  DynamicVMS<TDim>& rThis)
616 {
617  return rIStream;
618 }
619 
621 template< unsigned int TDim >
622 inline std::ostream& operator <<(std::ostream& rOStream,
623  const DynamicVMS<TDim>& rThis)
624 {
625  rThis.PrintInfo(rOStream);
626  rOStream << std::endl;
627  rThis.PrintData(rOStream);
628 
629  return rOStream;
630 }
632 
634 
635 
636 } // namespace Kratos.
637 
638 #endif // KRATOS_DYNAMIC_VMS__H_INCLUDED defined
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
A stabilized element for the incompressible Navier-Stokes equations.
Definition: dynamic_vms.h:106
virtual void UpdateSubscale(const ProcessInfo &rCurrentProcessInfo)
Calculate the value of the subscale velocity for next iteration.
Definition: dynamic_vms.cpp:617
void GetFirstDerivativesVector(Vector &rValues, int Step=0) const override
Returns VELOCITY_X, VELOCITY_Y, (VELOCITY_Z,) PRESSURE for each node.
Definition: dynamic_vms.cpp:332
void GetSecondDerivativesVector(Vector &rValues, int Step=0) const override
Returns ACCELERATION_X, ACCELERATION_Y, (ACCELERATION_Z,) 0 for each node.
Definition: dynamic_vms.cpp:355
virtual void CalculateASGSVelocityContribution(MatrixType &rDampingMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: dynamic_vms.cpp:799
virtual void OSSMomentumResidual(array_1d< double, 3 > &rResult, const double Density, const array_1d< double, 3 > &rConvVel, const ShapeFunctionsType &rN)
Definition: dynamic_vms.cpp:1277
Kratos::Vector ShapeFunctionsType
Type for shape function values container.
Definition: dynamic_vms.h:115
virtual void EvaluateViscosity(double &rViscosity, const ShapeFunctionsType &rN)
Evaluate kinematic viscosity at the given area coordinates. This function is intended to be used for ...
Definition: dynamic_vms.cpp:1198
void EvaluateVorticity(array_1d< double, 3 > &rVorticity, const ShapeDerivativesType &rDN_DX)
void Initialize(const ProcessInfo &rCurrentProcessInfo) override
Initialize containters for subscales on integration points.
Definition: dynamic_vms.cpp:90
std::string Info() const override
Turn back information as a string.
Definition: dynamic_vms.cpp:552
virtual double TauTime(const double Density, const double Viscosity, const double ConvVel, const double Dt)
Calculate 1 / ( rho/dt + 1/TauOne )
Definition: dynamic_vms.cpp:1231
void Calculate(const Variable< array_1d< double, 3 > > &rVariable, array_1d< double, 3 > &rOutput, const ProcessInfo &rCurrentProcessInfo) override
Definition: dynamic_vms.cpp:378
Kratos::Matrix ShapeDerivativesType
Type for shape function derivatives container.
Definition: dynamic_vms.h:118
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(DynamicVMS)
Pointer definition of DynamicVMS.
virtual void ConvectionOperator(Vector &rResult, const array_1d< double, 3 > &rConvVel)
Definition: dynamic_vms.cpp:1204
virtual void CalculateGeometryData()
Initialize shape functions derivatives and calculate the determinant of the element's Jacobian.
Definition: dynamic_vms.cpp:577
virtual void ASGSMomentumResidual(array_1d< double, 3 > &rResult, const double Density, const array_1d< double, 3 > &rConvVel, const ShapeFunctionsType &rN)
Definition: dynamic_vms.cpp:1249
~DynamicVMS() override
Destructor.
Definition: dynamic_vms.cpp:72
virtual void CalculateOSSVelocityContribution(MatrixType &rDampingMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: dynamic_vms.cpp:940
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: dynamic_vms.cpp:566
virtual double InvTauTime(const double Density, const double Viscosity, const double ConvVel, const double Dt)
Calculate 1 / TauTime.
Definition: dynamic_vms.cpp:1240
void EvaluateInPoint(TValueType &rValue, const Kratos::Variable< TValueType > &rVariable, const ShapeFunctionsType &rN)
Definition: dynamic_vms.h:391
virtual void DenseSystemSolve(const Matrix &rA, const Vector &rB, Vector &rX)
void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const override
Provides the global indices for each one of this element's local rows.
virtual void MassResidual(double &rResult)
Definition: dynamic_vms.cpp:1305
GeometryData::IntegrationMethod GetIntegrationMethod() const override
Accessor to the integration method.
Definition: dynamic_vms.cpp:545
virtual void LumpedMassMatrix(MatrixType &rMassMatrix)
Definition: dynamic_vms.cpp:1089
void EvaluateInPoint(TValueType &rValue, const Kratos::Variable< TValueType > &rVariable, const ShapeFunctionsType &rN, const unsigned int Step)
Definition: dynamic_vms.h:404
void CalculateLocalVelocityContribution(MatrixType &rDampingMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Computes the local contribution associated to 'new' velocity and pressure values.
Definition: dynamic_vms.cpp:149
void InitializeSolutionStep(const ProcessInfo &rCurrentProcessInfo) override
Prepare the element for a new solution step. Update the values on the subscales and evaluate elementa...
Definition: dynamic_vms.cpp:94
void InitializeNonLinearIteration(const ProcessInfo &rCurrentProcessInfo) override
Calculate a new value for the velocity subscale.
Definition: dynamic_vms.cpp:104
void CalculateOnIntegrationPoints(const Variable< Vector > &rVariable, std::vector< Vector > &rValues, const ProcessInfo &rCurrentProcessInfo) override
Empty implementation of unused CalculateOnIntegrationPoints overloads to avoid compilation warning.
Definition: dynamic_vms.h:301
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: dynamic_vms.cpp:560
void CalculateMassMatrix(MatrixType &rMassMatrix, const ProcessInfo &rCurrentProcessInfo) override
Definition: dynamic_vms.cpp:161
Element::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
Create a new element of this type.
Definition: dynamic_vms.cpp:77
void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: dynamic_vms.cpp:135
void CalculateOnIntegrationPoints(const Variable< Matrix > &rVariable, std::vector< Matrix > &rValues, const ProcessInfo &rCurrentProcessInfo) override
Empty implementation of unused CalculateOnIntegrationPoints overloads to avoid compilation warning.
Definition: dynamic_vms.h:307
void CalculateOnIntegrationPoints(const Variable< array_1d< double, 3 > > &rVariable, std::vector< array_1d< double, 3 > > &rOutput, const ProcessInfo &rCurrentProcessInfo) override
Definition: dynamic_vms.cpp:493
virtual void AddViscousTerm(MatrixType &rDampingMatrix, const double Weight, const ShapeDerivativesType &rDN_DX)
virtual double TauTwo(const double Density, const double Viscosity, const double ConvVel)
Calculate pressure subscale stabilization parameter.
Definition: dynamic_vms.cpp:1225
void CalculateOnIntegrationPoints(const Variable< array_1d< double, 6 > > &rVariable, std::vector< array_1d< double, 6 > > &rValues, const ProcessInfo &rCurrentProcessInfo) override
Empty implementation of unused CalculateOnIntegrationPoints overloads to avoid compilation warning.
Definition: dynamic_vms.h:295
virtual void LinearUpdateSubscale(const ProcessInfo &rCurrentProcessInfo)
Definition: dynamic_vms.cpp:747
void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: dynamic_vms.cpp:114
void SetValuesOnIntegrationPoints(const Variable< double > &rVariable, const std::vector< double > &rValues, const ProcessInfo &rCurrentProcessInfo) override
Definition: dynamic_vms.cpp:514
virtual void EvaluateConvVelocity(array_1d< double, 3 > &rConvection, const ShapeFunctionsType &rN)
Definition: dynamic_vms.cpp:1173
virtual void ConsistentMassMatrix(MatrixType &rMassMatrix)
Definition: dynamic_vms.cpp:1126
virtual double TauOne(const double Density, const double Viscosity, const double ConvVel)
Calculate velocity subscale stabilization parameter.
Definition: dynamic_vms.cpp:1216
void GetDofList(DofsVectorType &rElementalDofList, const ProcessInfo &rCurrentProcessInfo) const override
Returns a list of the element's Dofs.
Base class for all Elements.
Definition: element.h:60
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
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
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
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
Dt
Definition: face_heat.py:78
def load(f)
Definition: ode_solve.py:307
integer i
Definition: TensorModule.f:17