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.
alternative_qs_vms_dem_coupled.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: Joaquin Gonzalez-Usua
11 //
12 
13 #ifndef KRATOS_ALTERNATIVE_QS_VMS_DEM_COUPLED_H
14 #define KRATOS_ALTERNATIVE_QS_VMS_DEM_COUPLED_H
15 
16 // System includes
17 #include <string>
18 #include <iostream>
19 
20 // External includes
21 
22 // Project includes
23 #include "includes/checks.h"
24 #include "includes/element.h"
25 #include "includes/serializer.h"
26 #include "geometries/geometry.h"
27 #include "includes/cfd_variables.h"
28 
29 // Application includes
31 #include "custom_elements/qs_vms.h"
33 
34 namespace Kratos
35 {
36 
39 
42 
46 
50 
54 
58 
59 template< class TElementData >
60 class AlternativeQSVMSDEMCoupled : public QSVMS<TElementData>
61 {
62 public:
65 
68 
70  typedef Node NodeType;
71 
74 
77 
79  typedef Vector VectorType;
80 
82  typedef Matrix MatrixType;
83 
84  typedef std::size_t IndexType;
85 
86  typedef std::size_t SizeType;
87 
88  typedef std::vector<std::size_t> EquationIdVectorType;
89 
90  typedef std::vector< Dof<double>::Pointer > DofsVectorType;
91 
93 
96 
99 
102 
104 
105  constexpr static unsigned int Dim = QSVMS<TElementData>::Dim;
106  constexpr static unsigned int NumNodes = QSVMS<TElementData>::NumNodes;
107  constexpr static unsigned int BlockSize = QSVMS<TElementData>::BlockSize;
108  constexpr static unsigned int LocalSize = QSVMS<TElementData>::LocalSize;
109  constexpr static unsigned int StrainSize = QSVMS<TElementData>::StrainSize;
110 
114 
115  //Constructors.
116 
118 
122 
124 
128  AlternativeQSVMSDEMCoupled(IndexType NewId, const NodesArrayType& ThisNodes);
129 
131 
135  AlternativeQSVMSDEMCoupled(IndexType NewId, GeometryType::Pointer pGeometry);
136 
138 
143  AlternativeQSVMSDEMCoupled(IndexType NewId, GeometryType::Pointer pGeometry, Properties::Pointer pProperties);
144 
147 
151 
152 
156 
157 
159 
166  Element::Pointer Create(
167  IndexType NewId,
168  NodesArrayType const& ThisNodes,
169  Properties::Pointer pProperties) const override;
170 
172 
179  Element::Pointer Create(
180  IndexType NewId,
181  GeometryType::Pointer pGeom,
182  Properties::Pointer pProperties) const override;
183 
184 
185  void Initialize(const ProcessInfo& rCurrentProcessInfo) override;
186 
188 
192 
196 
197  int Check(const ProcessInfo &rCurrentProcessInfo) const override;
198 
202 
204  std::string Info() const override;
205 
206 
208  void PrintInfo(std::ostream& rOStream) const override;
209 
210 
214 
215 
217 
218 protected:
219 
222 
223 
227 
230  // Velocity subscale history, stored at integration points
236 
237 
241 
242  // Protected interface of FluidElement ////////////////////////////////////
243 
245 
247  const TElementData& rData,
248  const array_1d<double,3> &rConvectionVelocity,
249  array_1d<double,3>& rResidual) const override;
250 
251  void MomentumProjTerm(
252  const TElementData& rData,
253  const array_1d<double,3>& rConvectionVelocity,
254  array_1d<double,3> &rMomentumRHS) const override;
255 
257  TElementData& rData,
258  MatrixType &rMassMatrix);
259 
261  TElementData& rData,
262  BoundedMatrix<double,NumNodes*(Dim+1),NumNodes*(Dim+1)>& rLHS,
263  VectorType& rLocalRHS);
264 
265  void AddViscousTerm(
266  const TElementData& rData,
268  VectorType& rRHS) override;
269 
271  void CalculateTau(
272  const TElementData& rData,
273  const array_1d<double,3> &Velocity,
275  double &TauTwo) const;
276 
278  const ProcessInfo &rCurrentProcessInfo) override;
279 
281  TElementData& rData,
282  unsigned int IntegrationPointIndex,
283  double Weight,
284  const typename TElementData::MatrixRowType& rN,
285  const typename TElementData::ShapeDerivativesType& rDN_DX,
286  const typename TElementData::ShapeFunctionsSecondDerivativesType& rDDN_DDX) const;
287 
288  void AddVelocitySystem(
289  TElementData& rData,
290  MatrixType &rLocalLHS,
291  VectorType &rLocalRHS) override;
292 
293  void CalculateMassMatrix(MatrixType& rMassMatrix,
294  const ProcessInfo& rCurrentProcessInfo) override;
295 
297  VectorType& rRightHandSideVector,
298  const ProcessInfo& rCurrentProcessInfo) override;
299 
300  void AddMassLHS(
301  TElementData& rData,
302  MatrixType& rMassMatrix) override;
303 
305  const TElementData& rData);
306 
307  void MassProjTerm(
308  const TElementData& rData,
309  double &rMassRHS) const override;
310 
311  void SubscaleVelocity(
312  const TElementData& rData,
313  array_1d<double,3> &rVelocitySubscale) const override;
314 
315  void SubscalePressure(
316  const TElementData& rData,
317  double &rPressureSubscale) const override;
318 
319  void Calculate(
320  const Variable<array_1d<double, 3>>& rVariable,
321  array_1d<double, 3>& rOutput, const ProcessInfo& rCurrentProcessInfo) override;
322 
324  const Variable<array_1d<double, 3>>& rVariable,
325  std::vector<array_1d<double, 3>>& rOutput,
326  const ProcessInfo& rCurrentProcessInfo) override;
327 
329  const Variable<double>& rVariable,
330  std::vector<double>& rOutput,
331  const ProcessInfo& rCurrentProcessInfo) override;
332 
334  Variable<Matrix> const& rVariable,
335  std::vector<Matrix>& rValues,
336  ProcessInfo const& rCurrentProcessInfo) override;
337 
338  void InitializeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo) override;
339 
340  void FinalizeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo) override;
341 
342  void UpdateSubscaleVelocity(const TElementData& rData);
343 
344  array_1d<double,3> FullConvectiveVelocity(const TElementData& rData) const;
345 
349 
350 
354 
355 
359 
361 
362 private:
363 
366 
370 
374 
375  friend class Serializer;
376 
377  void save(Serializer& rSerializer) const override;
378 
379  void load(Serializer& rSerializer) override;
380 
384 
385 
389 
390 
394 
395 
399 
400 
404 
407 
411 
412 
413 }; // Class AlternativeQSVMSDEMCoupled
414 
416 
419 
420 
424 
425 
427 template< class TElementData >
428 inline std::istream& operator >>(std::istream& rIStream,
430 {
431  return rIStream;
432 }
433 
435 template< class TElementData >
436 inline std::ostream& operator <<(std::ostream& rOStream,
438 {
439  rThis.PrintInfo(rOStream);
440  rOStream << std::endl;
441  rThis.PrintData(rOStream);
442 
443  return rOStream;
444 }
446 
448 
449 } // namespace Kratos.
450 
451 #endif // KRATOS_ALTERNATIVE_QS_VMS_DEM_COUPLED_H
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
Definition: alternative_qs_vms_dem_coupled.h:61
void CalculateProjections(const ProcessInfo &rCurrentProcessInfo) override
Definition: alternative_qs_vms_dem_coupled.cpp:1151
Geometry< NodeType > GeometryType
Geometry type (using with given NodeType)
Definition: alternative_qs_vms_dem_coupled.h:73
PointerVectorSet< Dof< double >, IndexedObject > DofsArrayType
Definition: alternative_qs_vms_dem_coupled.h:92
void Initialize(const ProcessInfo &rCurrentProcessInfo) override
Definition: alternative_qs_vms_dem_coupled.cpp:70
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(AlternativeQSVMSDEMCoupled)
Pointer definition of AlternativeQSVMSDEMCoupled.
void CalculateResistanceTensor(const TElementData &rData)
Definition: alternative_qs_vms_dem_coupled.cpp:1041
void SubscaleVelocity(const TElementData &rData, array_1d< double, 3 > &rVelocitySubscale) const override
Definition: alternative_qs_vms_dem_coupled.cpp:1210
constexpr static unsigned int BlockSize
Definition: alternative_qs_vms_dem_coupled.h:107
void MassProjTerm(const TElementData &rData, double &rMassRHS) const override
Definition: alternative_qs_vms_dem_coupled.cpp:1082
constexpr static unsigned int Dim
Definition: alternative_qs_vms_dem_coupled.h:105
std::vector< Dof< double >::Pointer > DofsVectorType
Definition: alternative_qs_vms_dem_coupled.h:90
std::string Info() const override
Turn back information as a string.
Definition: alternative_qs_vms_dem_coupled.cpp:338
~AlternativeQSVMSDEMCoupled()
Destructor.
Definition: alternative_qs_vms_dem_coupled.cpp:54
void AddMassStabilization(TElementData &rData, MatrixType &rMassMatrix)
Definition: alternative_qs_vms_dem_coupled.cpp:508
std::vector< std::size_t > EquationIdVectorType
Definition: alternative_qs_vms_dem_coupled.h:88
GeometryData::IntegrationMethod GetIntegrationMethod() const override
Definition: alternative_qs_vms_dem_coupled.cpp:306
constexpr static unsigned int LocalSize
Definition: alternative_qs_vms_dem_coupled.h:108
DenseVector< array_1d< double, Dim > > mPredictedSubscaleVelocity
Definition: alternative_qs_vms_dem_coupled.h:231
void UpdateIntegrationPointDataSecondDerivatives(TElementData &rData, unsigned int IntegrationPointIndex, double Weight, const typename TElementData::MatrixRowType &rN, const typename TElementData::ShapeDerivativesType &rDN_DX, const typename TElementData::ShapeFunctionsSecondDerivativesType &rDDN_DDX) const
Definition: alternative_qs_vms_dem_coupled.cpp:398
constexpr static unsigned int NumNodes
Definition: alternative_qs_vms_dem_coupled.h:106
Element::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, Properties::Pointer pProperties) const override
Create a new element of this type.
Definition: alternative_qs_vms_dem_coupled.cpp:58
Kratos::Vector ShapeFunctionsType
Type for shape function values container.
Definition: alternative_qs_vms_dem_coupled.h:95
void UpdateSubscaleVelocity(const TElementData &rData)
Definition: alternative_qs_vms_dem_coupled.cpp:1270
AlternativeQSVMSDEMCoupled(IndexType NewId=0)
Default constuctor.
Definition: alternative_qs_vms_dem_coupled.cpp:30
std::size_t SizeType
Definition: alternative_qs_vms_dem_coupled.h:86
array_1d< double, 3 > FullConvectiveVelocity(const TElementData &rData) const
Definition: alternative_qs_vms_dem_coupled.cpp:1255
Vector VectorType
Vector type for local contributions to the linear system.
Definition: alternative_qs_vms_dem_coupled.h:79
DenseVector< array_1d< double, Dim > > mPreviousVelocity
Definition: alternative_qs_vms_dem_coupled.h:232
void CalculateOnIntegrationPoints(const Variable< array_1d< double, 3 >> &rVariable, std::vector< array_1d< double, 3 >> &rOutput, const ProcessInfo &rCurrentProcessInfo) override
Definition: alternative_qs_vms_dem_coupled.cpp:181
void AddViscousTerm(const TElementData &rData, BoundedMatrix< double, LocalSize, LocalSize > &rLHS, VectorType &rRHS) override
Definition: alternative_qs_vms_dem_coupled.cpp:667
void CalculateMassMatrix(MatrixType &rMassMatrix, const ProcessInfo &rCurrentProcessInfo) override
Definition: alternative_qs_vms_dem_coupled.cpp:963
void Calculate(const Variable< array_1d< double, 3 >> &rVariable, array_1d< double, 3 > &rOutput, const ProcessInfo &rCurrentProcessInfo) override
Definition: alternative_qs_vms_dem_coupled.cpp:224
Kratos::Matrix ShapeFunctionDerivativesType
Type for a matrix containing the shape function gradients.
Definition: alternative_qs_vms_dem_coupled.h:98
Matrix MatrixType
Matrix type for local contributions to the linear system.
Definition: alternative_qs_vms_dem_coupled.h:82
void AddReactionStabilization(TElementData &rData, BoundedMatrix< double, NumNodes *(Dim+1), NumNodes *(Dim+1)> &rLHS, VectorType &rLocalRHS)
Definition: alternative_qs_vms_dem_coupled.cpp:572
GeometryType::ShapeFunctionsGradientsType ShapeFunctionDerivativesArrayType
Type for an array of shape function gradient matrices.
Definition: alternative_qs_vms_dem_coupled.h:101
Geometry< NodeType >::PointsArrayType NodesArrayType
Definition of nodes container type, redefined from GeometryType.
Definition: alternative_qs_vms_dem_coupled.h:76
void AlgebraicMomentumResidual(const TElementData &rData, const array_1d< double, 3 > &rConvectionVelocity, array_1d< double, 3 > &rResidual) const override
Determine the shape second derivative in the gauss point.
Definition: alternative_qs_vms_dem_coupled.cpp:411
void MomentumProjTerm(const TElementData &rData, const array_1d< double, 3 > &rConvectionVelocity, array_1d< double, 3 > &rMomentumRHS) const override
Definition: alternative_qs_vms_dem_coupled.cpp:462
void AddMassLHS(TElementData &rData, MatrixType &rMassMatrix) override
Definition: alternative_qs_vms_dem_coupled.cpp:1051
void SubscalePressure(const TElementData &rData, double &rPressureSubscale) const override
Definition: alternative_qs_vms_dem_coupled.cpp:1233
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: alternative_qs_vms_dem_coupled.cpp:347
Node NodeType
Node type (default is: Node)
Definition: alternative_qs_vms_dem_coupled.h:70
void CalculateLocalVelocityContribution(MatrixType &rDampMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: alternative_qs_vms_dem_coupled.cpp:999
void CalculateTau(const TElementData &rData, const array_1d< double, 3 > &Velocity, BoundedMatrix< double, Dim, Dim > &TauOne, double &TauTwo) const
Definition: alternative_qs_vms_dem_coupled.cpp:1102
GeometryType::ShapeFunctionsSecondDerivativesType ShapeFunctionsSecondDerivativesType
Definition: alternative_qs_vms_dem_coupled.h:103
void AddVelocitySystem(TElementData &rData, MatrixType &rLocalLHS, VectorType &rLocalRHS) override
Definition: alternative_qs_vms_dem_coupled.cpp:688
int mInterpolationOrder
Definition: alternative_qs_vms_dem_coupled.h:228
constexpr static unsigned int StrainSize
Definition: alternative_qs_vms_dem_coupled.h:109
DenseVector< BoundedMatrix< double, Dim, Dim > > mViscousResistanceTensor
Definition: alternative_qs_vms_dem_coupled.h:229
void InitializeNonLinearIteration(const ProcessInfo &rCurrentProcessInfo) override
Definition: alternative_qs_vms_dem_coupled.cpp:353
std::size_t IndexType
Definition: alternative_qs_vms_dem_coupled.h:84
void FinalizeNonLinearIteration(const ProcessInfo &rCurrentProcessInfo) override
Definition: alternative_qs_vms_dem_coupled.cpp:376
int Check(const ProcessInfo &rCurrentProcessInfo) const override
Definition: alternative_qs_vms_dem_coupled.cpp:318
std::size_t IndexType
Definition: flags.h:74
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: fluid_element.hpp:584
IntegrationMethod
Definition: geometry_data.h:76
Geometry base class.
Definition: geometry.h:71
This object defines an indexed object.
Definition: indexed_object.h:54
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
A sorted associative container similar to an STL set, but uses a vector to store pointers to its data...
Definition: pointer_vector_set.h:72
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
Definition: qs_vms.h:52
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
def load(f)
Definition: ode_solve.py:307