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.
d_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_D_VMS_DEM_COUPLED_H
14 #define KRATOS_D_VMS_DEM_COUPLED_H
15 
16 #include "includes/define.h"
17 #include "includes/element.h"
18 #include "includes/serializer.h"
19 #include "geometries/geometry.h"
20 
21 #include "includes/cfd_variables.h"
23 #include "custom_elements/d_vms.h"
25 
26 namespace Kratos
27 {
28 
31 
34 
38 
42 
46 
50 
51 template< class TElementData >
52 class DVMSDEMCoupled : public DVMS<TElementData>
53 {
54 public:
57 
60 
62  typedef Node NodeType;
63 
66 
69 
71  typedef Vector VectorType;
72 
74  typedef Matrix MatrixType;
75 
76  typedef std::size_t IndexType;
77 
78  typedef std::size_t SizeType;
79 
80  typedef std::vector<std::size_t> EquationIdVectorType;
81 
82  typedef std::vector< Dof<double>::Pointer > DofsVectorType;
83 
85 
88 
91 
94 
96 
97  constexpr static unsigned int Dim = DVMS<TElementData>::Dim;
98  constexpr static unsigned int NumNodes = DVMS<TElementData>::NumNodes;
99  constexpr static unsigned int BlockSize = DVMS<TElementData>::BlockSize;
100  constexpr static unsigned int LocalSize = DVMS<TElementData>::LocalSize;
101  constexpr static unsigned int StrainSize = DVMS<TElementData>::StrainSize;
102 
106 
107  //Constructors.
108 
110 
113  DVMSDEMCoupled(IndexType NewId = 0);
114 
116 
120  DVMSDEMCoupled(IndexType NewId, const NodesArrayType& ThisNodes);
121 
123 
127  DVMSDEMCoupled(IndexType NewId, GeometryType::Pointer pGeometry);
128 
130 
135  DVMSDEMCoupled(IndexType NewId, GeometryType::Pointer pGeometry, Properties::Pointer pProperties);
136 
138  ~DVMSDEMCoupled();
139 
143 
144 
148 
149 
151 
158  Element::Pointer Create(
159  IndexType NewId,
160  NodesArrayType const& ThisNodes,
161  Properties::Pointer pProperties) const override;
162 
164 
171  Element::Pointer Create(
172  IndexType NewId,
173  GeometryType::Pointer pGeom,
174  Properties::Pointer pProperties) const override;
175 
176  void Initialize(const ProcessInfo& rCurrentProcessInfo) override;
177 
179 
181  const Variable<array_1d<double, 3>>& rVariable,
182  std::vector<array_1d<double, 3>>& rOutput,
183  const ProcessInfo& rCurrentProcessInfo) override;
184 
186  const Variable<double>& rVariable,
187  std::vector<double>& rOutput,
188  const ProcessInfo& rCurrentProcessInfo) override;
189 
191  Variable<Matrix> const& rVariable,
192  std::vector<Matrix>& rValues,
193  ProcessInfo const& rCurrentProcessInfo) override;
194 
195  void FinalizeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
196 
197  void InitializeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo) override;
198 
199  void FinalizeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo) override;
200 
204 
208 
210  std::string Info() const override;
211 
212 
214  void PrintInfo(std::ostream& rOStream) const override;
215 
216 
220 
221 
223 
224 protected:
225 
228 
234  // Velocity subscale history, stored at integration points
238 
242 
243 
247 
248  // Protected interface of FluidElement ////////////////////////////////////
249 
251  const TElementData& rData,
252  const array_1d<double,3> &rConvectionVelocity,
253  array_1d<double,3>& rResidual) const override;
254 
255  void MomentumProjTerm(
256  const TElementData& rData,
257  const array_1d<double,3>& rConvectionVelocity,
258  array_1d<double,3> &rMomentumRHS) const override;
259 
260  void AddVelocitySystem(
261  TElementData& rData,
262  MatrixType& rLocalLHS,
263  VectorType& rLocalRHS) override;
264 
266  TElementData& rData,
267  BoundedMatrix<double,NumNodes*(Dim+1),NumNodes*(Dim+1)>& rLHS,
268  VectorType& rLocalRHS);
269 
271  TElementData& rData,
272  MatrixType& rMassMatrix) override;
273 
274  void CalculateProjections(const ProcessInfo &rCurrentProcessInfo) override;
275 
277  const TElementData& rData,
278  const array_1d<double,3> &Velocity,
280  double &TauTwo) const;
281 
283  TElementData& rData,
284  unsigned int IntegrationPointIndex,
285  double Weight,
286  const typename TElementData::MatrixRowType& rN,
287  const typename TElementData::ShapeDerivativesType& rDN_DX,
288  const typename TElementData::ShapeFunctionsSecondDerivativesType& rDDN_DDX) const;
289 
290  void SubscaleVelocity(
291  const TElementData& rData,
292  array_1d<double,3>& rVelocitySubscale) const override;
293 
294  void SubscalePressure(
295  const TElementData& rData,
296  double& rPressureSubscale) const override;
297 
298  void CalculateMassMatrix(MatrixType& rMassMatrix,
299  const ProcessInfo& rCurrentProcessInfo) override;
300 
302  VectorType& rRightHandSideVector,
303  const ProcessInfo& rCurrentProcessInfo) override;
304 
306  const TElementData& rData) const override;
307 
309  const TElementData& rData) override;
310 
312  const TElementData& rData);
313 
315  const TElementData& rData);
316 
317  void AddMassLHS(
318  TElementData& rData,
319  MatrixType& rMassMatrix) override;
320 
321  void MassProjTerm(
322  const TElementData& rData,
323  double& rMassRHS) const override;
324 
328 
329 
333 
334 
338 
340 
341 private:
342 
345 
349 
350  friend class Serializer;
351 
352  void save(Serializer& rSerializer) const override;
353 
354  void load(Serializer& rSerializer) override;
355 
359 
360 
364 
368 
369 
373 
374 
378 
380  DVMSDEMCoupled& operator=(DVMSDEMCoupled const& rOther);
381 
383  DVMSDEMCoupled(DVMSDEMCoupled const& rOther);
384 
386 
387 
388 }; // Class DVMSDEMCoupled
389 
391 
394 
395 
399 
400 
402 template< class TElementData >
403 inline std::istream& operator >>(std::istream& rIStream,
405 {
406  return rIStream;
407 }
408 
410 template< class TElementData >
411 inline std::ostream& operator <<(std::ostream& rOStream,
412  const DVMSDEMCoupled<TElementData>& rThis)
413 {
414  rThis.PrintInfo(rOStream);
415  rOStream << std::endl;
416  rThis.PrintData(rOStream);
417 
418  return rOStream;
419 }
421 
423 
424 } // namespace Kratos.
425 
426 #endif // KRATOS_D_VMS_DEM_COUPLED_H
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
Definition: d_vms_dem_coupled.h:53
constexpr static unsigned int Dim
Definition: d_vms_dem_coupled.h:97
Node NodeType
Node type (default is: Node)
Definition: d_vms_dem_coupled.h:62
void Initialize(const ProcessInfo &rCurrentProcessInfo) override
Set up the element.
Definition: d_vms_dem_coupled.cpp:207
GeometryType::ShapeFunctionsGradientsType ShapeFunctionDerivativesArrayType
Type for an array of shape function gradient matrices.
Definition: d_vms_dem_coupled.h:93
std::string Info() const override
Turn back information as a string.
Definition: d_vms_dem_coupled.cpp:361
DenseVector< array_1d< double, Dim > > mOldSubscaleVelocity
Definition: d_vms_dem_coupled.h:236
void UpdateSubscaleVelocity(const TElementData &rData)
Definition: d_vms_dem_coupled.cpp:1124
constexpr static unsigned int LocalSize
Definition: d_vms_dem_coupled.h:100
void MassProjTerm(const TElementData &rData, double &rMassRHS) const override
Definition: d_vms_dem_coupled.cpp:869
void UpdateSubscaleVelocityPrediction(const TElementData &rData) override
Definition: d_vms_dem_coupled.cpp:1163
PointerVectorSet< Dof< double >, IndexedObject > DofsArrayType
Definition: d_vms_dem_coupled.h:84
~DVMSDEMCoupled()
Destructor.
Definition: d_vms_dem_coupled.cpp:68
constexpr static unsigned int StrainSize
Definition: d_vms_dem_coupled.h:101
std::vector< Dof< double >::Pointer > DofsVectorType
Definition: d_vms_dem_coupled.h:82
void FinalizeNonLinearIteration(const ProcessInfo &rCurrentProcessInfo) override
Definition: d_vms_dem_coupled.cpp:297
Matrix MatrixType
Matrix type for local contributions to the linear system.
Definition: d_vms_dem_coupled.h:74
array_1d< double, 3 > FullConvectiveVelocity(const TElementData &rData) const override
Definition: d_vms_dem_coupled.cpp:1109
std::size_t IndexType
Definition: d_vms_dem_coupled.h:76
DenseVector< array_1d< double, Dim > > mPreviousVelocity
Definition: d_vms_dem_coupled.h:237
GeometryType::ShapeFunctionsSecondDerivativesType ShapeFunctionsSecondDerivativesType
Definition: d_vms_dem_coupled.h:95
Kratos::Vector ShapeFunctionsType
Type for shape function values container.
Definition: d_vms_dem_coupled.h:87
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(DVMSDEMCoupled)
Pointer definition of DVMSDEMCoupled.
Vector VectorType
Vector type for local contributions to the linear system.
Definition: d_vms_dem_coupled.h:71
std::vector< std::size_t > EquationIdVectorType
Definition: d_vms_dem_coupled.h:80
void FinalizeSolutionStep(const ProcessInfo &rCurrentProcessInfo) override
Update the values of tracked small scale quantities.
Definition: d_vms_dem_coupled.cpp:332
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: d_vms_dem_coupled.cpp:370
GeometryData::IntegrationMethod GetIntegrationMethod() const override
Definition: d_vms_dem_coupled.cpp:266
constexpr static unsigned int BlockSize
Definition: d_vms_dem_coupled.h:99
std::size_t SizeType
Definition: d_vms_dem_coupled.h:78
void CalculateMassMatrix(MatrixType &rMassMatrix, const ProcessInfo &rCurrentProcessInfo) override
Definition: d_vms_dem_coupled.cpp:749
void CalculateStabilizationParameters(const TElementData &rData, const array_1d< double, 3 > &Velocity, BoundedMatrix< double, Dim, Dim > &TauOne, double &TauTwo) const
Definition: d_vms_dem_coupled.cpp:1013
void MomentumProjTerm(const TElementData &rData, const array_1d< double, 3 > &rConvectionVelocity, array_1d< double, 3 > &rMomentumRHS) const override
Definition: d_vms_dem_coupled.cpp:424
DenseVector< array_1d< double, Dim > > mPredictedSubscaleVelocity
Definition: d_vms_dem_coupled.h:235
void AlgebraicMomentumResidual(const TElementData &rData, const array_1d< double, 3 > &rConvectionVelocity, array_1d< double, 3 > &rResidual) const override
Definition: d_vms_dem_coupled.cpp:382
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: d_vms_dem_coupled.cpp:319
Geometry< NodeType > GeometryType
Geometry type (using with given NodeType)
Definition: d_vms_dem_coupled.h:65
Kratos::Matrix ShapeFunctionDerivativesType
Type for a matrix containing the shape function gradients.
Definition: d_vms_dem_coupled.h:90
void SubscalePressure(const TElementData &rData, double &rPressureSubscale) const override
Definition: d_vms_dem_coupled.cpp:1088
void AddVelocitySystem(TElementData &rData, MatrixType &rLocalLHS, VectorType &rLocalRHS) override
Definition: d_vms_dem_coupled.cpp:548
constexpr static unsigned int NumNodes
Definition: d_vms_dem_coupled.h:98
void SubscaleVelocity(const TElementData &rData, array_1d< double, 3 > &rVelocitySubscale) const override
Definition: d_vms_dem_coupled.cpp:1058
int mInterpolationOrder
Definition: d_vms_dem_coupled.h:232
Geometry< NodeType >::PointsArrayType NodesArrayType
Definition of nodes container type, redefined from GeometryType.
Definition: d_vms_dem_coupled.h:68
void AddMassStabilization(TElementData &rData, MatrixType &rMassMatrix) override
Definition: d_vms_dem_coupled.cpp:893
DVMSDEMCoupled(IndexType NewId=0)
Default constuctor.
Definition: d_vms_dem_coupled.cpp:32
DenseVector< BoundedMatrix< double, Dim, Dim > > mViscousResistanceTensor
Definition: d_vms_dem_coupled.h:233
void CalculateProjections(const ProcessInfo &rCurrentProcessInfo) override
Definition: d_vms_dem_coupled.cpp:953
void InitializeNonLinearIteration(const ProcessInfo &rCurrentProcessInfo) override
Predict the value of the small scale velocity for the current iteration.
Definition: d_vms_dem_coupled.cpp:275
void AddMassLHS(TElementData &rData, MatrixType &rMassMatrix) override
Definition: d_vms_dem_coupled.cpp:838
void CalculateLocalVelocityContribution(MatrixType &rDampMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: d_vms_dem_coupled.cpp:785
void CalculateResistanceTensor(const TElementData &rData)
Definition: d_vms_dem_coupled.cpp:828
void AddReactionStabilization(TElementData &rData, BoundedMatrix< double, NumNodes *(Dim+1), NumNodes *(Dim+1)> &rLHS, VectorType &rLocalRHS)
Definition: d_vms_dem_coupled.cpp:459
void CalculateOnIntegrationPoints(const Variable< array_1d< double, 3 >> &rVariable, std::vector< array_1d< double, 3 >> &rOutput, const ProcessInfo &rCurrentProcessInfo) override
Definition: d_vms_dem_coupled.cpp:164
Element::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, Properties::Pointer pProperties) const override
Create a new element of this type.
Definition: d_vms_dem_coupled.cpp:75
Definition: d_vms.h:54
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
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