13 #ifndef KRATOS_QS_VMS_DEM_COUPLED_H
14 #define KRATOS_QS_VMS_DEM_COUPLED_H
59 template<
class TElementData >
169 Properties::Pointer pProperties)
const override;
181 GeometryType::Pointer pGeom,
182 Properties::Pointer pProperties)
const override;
205 std::vector<double>& rOutput,
210 std::vector<Matrix>& rValues,
221 const ProcessInfo& rCurrentProcessInfo)
const override;
242 std::string
Info()
const override;
246 void PrintInfo(std::ostream& rOStream)
const override;
283 const TElementData& rData,
288 const TElementData& rData,
303 const TElementData& rData,
306 double &TauTwo)
const;
313 unsigned int IntegrationPointIndex,
315 const typename TElementData::MatrixRowType& rN,
316 const typename TElementData::ShapeDerivativesType& rDN_DX,
317 const typename TElementData::ShapeFunctionsSecondDerivativesType& rDDN_DDX)
const;
332 const TElementData& rData);
339 const TElementData& rData,
340 double &rMassRHS)
const override;
343 const TElementData& rData,
347 const TElementData& rData,
348 double &rPressureSubscale)
const override;
381 void save(
Serializer& rSerializer)
const override;
431 template<
class TElementData >
439 template<
class TElementData >
444 rOStream << std::endl;
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
std::vector< std::size_t > EquationIdVectorType
Definition: element.h:98
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_dem_coupled.h:61
constexpr static unsigned int NumNodes
Definition: qs_vms_dem_coupled.h:106
constexpr static unsigned int Dim
Definition: qs_vms_dem_coupled.h:105
constexpr static unsigned int BlockSize
Definition: qs_vms_dem_coupled.h:107
DenseVector< array_1d< double, Dim > > mPreviousVelocity
Definition: qs_vms_dem_coupled.h:269
void CalculateOnIntegrationPoints(const Variable< array_1d< double, 3 >> &rVariable, std::vector< array_1d< double, 3 >> &rOutput, const ProcessInfo &rCurrentProcessInfo) override
Definition: qs_vms_dem_coupled.cpp:136
void CalculateLocalVelocityContribution(MatrixType &rDampMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
CalculateLocalVelocityContribution Calculate the local contribution in terms of velocity and pressure...
Definition: qs_vms_dem_coupled.cpp:807
void AddReactionStabilization(TElementData &rData, BoundedMatrix< double, NumNodes *(Dim+1), NumNodes *(Dim+1)> &rLHS, VectorType &rLocalRHS)
Definition: qs_vms_dem_coupled.cpp:452
void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const override
Definition: qs_vms_dem_coupled.cpp:280
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(QSVMSDEMCoupled)
Pointer definition of QSVMSDEMCoupled.
constexpr static unsigned int StrainSize
Definition: qs_vms_dem_coupled.h:109
void CalculateTau(const TElementData &rData, const array_1d< double, 3 > &Velocity, BoundedMatrix< double, Dim, Dim > &TauOne, double &TauTwo) const
Definition: qs_vms_dem_coupled.cpp:912
GeometryType::ShapeFunctionsGradientsType ShapeFunctionDerivativesArrayType
Type for an array of shape function gradient matrices.
Definition: qs_vms_dem_coupled.h:101
void FinalizeNonLinearIteration(const ProcessInfo &rCurrentProcessInfo) override
Definition: qs_vms_dem_coupled.cpp:326
GeometryData::IntegrationMethod GetIntegrationMethod() const override
GetIntegrationMethod Return the integration order to be used.
Definition: qs_vms_dem_coupled.cpp:227
Matrix MatrixType
Matrix type for local contributions to the linear system.
Definition: qs_vms_dem_coupled.h:82
Geometry< NodeType > GeometryType
Geometry type (using with given NodeType)
Definition: qs_vms_dem_coupled.h:73
void AddVelocitySystem(TElementData &rData, MatrixType &rLocalLHS, VectorType &rLocalRHS) override
Definition: qs_vms_dem_coupled.cpp:595
~QSVMSDEMCoupled() override
Destructor.
Definition: qs_vms_dem_coupled.cpp:54
void Calculate(const Variable< double > &rVariable, double &rOutput, const ProcessInfo &rCurrentProcessInfo) override
Definition: qs_vms_dem_coupled.cpp:258
void CalculateResistanceTensor(const TElementData &rData)
Definition: qs_vms_dem_coupled.cpp:849
Kratos::Matrix ShapeFunctionDerivativesType
Type for a matrix containing the shape function gradients.
Definition: qs_vms_dem_coupled.h:98
void AddMassStabilization(TElementData &rData, MatrixType &rMassMatrix)
Definition: qs_vms_dem_coupled.cpp:536
std::size_t SizeType
Definition: qs_vms_dem_coupled.h:86
Element::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, Properties::Pointer pProperties) const override
Create a new element of this type.
Definition: qs_vms_dem_coupled.cpp:58
void SubscalePressure(const TElementData &rData, double &rPressureSubscale) const override
Definition: qs_vms_dem_coupled.cpp:1036
std::string Info() const override
Turn back information as a string.
Definition: qs_vms_dem_coupled.cpp:288
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: qs_vms_dem_coupled.cpp:348
DenseVector< BoundedMatrix< double, Dim, Dim > > mViscousResistanceTensor
Definition: qs_vms_dem_coupled.h:266
int Check(const ProcessInfo &rCurrentProcessInfo) const override
Definition: qs_vms_dem_coupled.cpp:239
std::vector< Dof< double >::Pointer > DofsVectorType
Definition: qs_vms_dem_coupled.h:90
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: qs_vms_dem_coupled.cpp:297
DenseVector< array_1d< double, Dim > > mPredictedSubscaleVelocity
Definition: qs_vms_dem_coupled.h:268
int mInterpolationOrder
Definition: qs_vms_dem_coupled.h:265
void AlgebraicMomentumResidual(const TElementData &rData, const array_1d< double, 3 > &rConvectionVelocity, array_1d< double, 3 > &rResidual) const override
Definition: qs_vms_dem_coupled.cpp:373
void Initialize(const ProcessInfo &rCurrentProcessInfo) override
Set up the element for solution.
Definition: qs_vms_dem_coupled.cpp:180
Kratos::Vector ShapeFunctionsType
Type for shape function values container.
Definition: qs_vms_dem_coupled.h:95
constexpr static unsigned int LocalSize
Definition: qs_vms_dem_coupled.h:108
void MomentumProjTerm(const TElementData &rData, const array_1d< double, 3 > &rConvectionVelocity, array_1d< double, 3 > &rMomentumRHS) const override
Definition: qs_vms_dem_coupled.cpp:412
QSVMSDEMCoupled(IndexType NewId=0)
Default constuctor.
Definition: qs_vms_dem_coupled.cpp:30
std::size_t IndexType
Definition: qs_vms_dem_coupled.h:84
void AddMassLHS(TElementData &rData, MatrixType &rMassMatrix) override
Definition: qs_vms_dem_coupled.cpp:860
void InitializeNonLinearIteration(const ProcessInfo &rCurrentProcessInfo) override
Definition: qs_vms_dem_coupled.cpp:303
Vector VectorType
Vector type for local contributions to the linear system.
Definition: qs_vms_dem_coupled.h:79
void UpdateSubscaleVelocity(const TElementData &rData)
Definition: qs_vms_dem_coupled.cpp:1058
void MassProjTerm(const TElementData &rData, double &rMassRHS) const override
Definition: qs_vms_dem_coupled.cpp:890
void CalculateMassMatrix(MatrixType &rMassMatrix, const ProcessInfo &rCurrentProcessInfo) override
MassMatrix Calculate the local mass matrix.
Definition: qs_vms_dem_coupled.cpp:771
void CalculateProjections(const ProcessInfo &rCurrentProcessInfo) override
Definition: qs_vms_dem_coupled.cpp:956
void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
CalculateRightHandSide Return an empty matrix of appropriate size. This element does not have a local...
Definition: qs_vms_dem_coupled.cpp:361
PointerVectorSet< Dof< double >, IndexedObject > DofsArrayType
Definition: qs_vms_dem_coupled.h:92
Node NodeType
Node type (default is: Node)
Definition: qs_vms_dem_coupled.h:70
GeometryType::ShapeFunctionsSecondDerivativesType ShapeFunctionsSecondDerivativesType
Definition: qs_vms_dem_coupled.h:103
Geometry< NodeType >::PointsArrayType NodesArrayType
Definition of nodes container type, redefined from GeometryType.
Definition: qs_vms_dem_coupled.h:76
std::vector< std::size_t > EquationIdVectorType
Definition: qs_vms_dem_coupled.h:88
void SubscaleVelocity(const TElementData &rData, array_1d< double, 3 > &rVelocitySubscale) const override
Definition: qs_vms_dem_coupled.cpp:1015
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
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