13 #ifndef KRATOS_D_CONVECTION_DIFFUSION_EXPLICIT_H
14 #define KRATOS_D_CONVECTION_DIFFUSION_EXPLICIT_H
61 template<
unsigned int TDim,
unsigned int TNumNodes>
90 GeometryType::Pointer pGeometry);
94 GeometryType::Pointer pGeometry,
95 Properties::Pointer pProperties);
114 Properties::Pointer pProperties)
const override;
118 GeometryType::Pointer pGeom,
119 Properties::Pointer pProperties)
const override;
150 std::string
Info()
const override
152 return "DConvectionDiffusionExplicitElement #";
158 rOStream <<
Info() << this->
Id();
213 void save(
Serializer& rSerializer)
const override
236 void DCalculateRightHandSideInternal(
237 BoundedVector<double, TNumNodes>& rRightHandSideBoundedVector,
240 void DCalculateOrthogonalSubgridScaleRHSInternal(
241 BoundedVector<double, TNumNodes>& rRightHandSideVector,
244 void UpdateUnknownSubgridScaleGaussPoint(
280 template<
unsigned int TDim,
unsigned int TNumNodes = TDim + 1>
282 std::istream& rIStream,
289 template<
unsigned int TDim,
unsigned int TNumNodes = TDim + 1>
291 std::ostream& rOStream,
295 rOStream << std::endl;
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
This element solves the convection-diffusion equation, stabilized with algebraic subgrid scale or ort...
Definition: d_convection_diffusion_explicit.h:63
Vector VectorType
Definition: d_convection_diffusion_explicit.h:73
virtual ~DConvectionDiffusionExplicit()
Default constuctor.
Definition: d_convection_diffusion_explicit.cpp:39
Node NodeType
Definition: d_convection_diffusion_explicit.h:70
Geometry< NodeType >::PointsArrayType NodesArrayType
Definition: d_convection_diffusion_explicit.h:72
BaseType::ElementData ElementData
Definition: d_convection_diffusion_explicit.h:69
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: d_convection_diffusion_explicit.h:156
Geometry< NodeType > GeometryType
Definition: d_convection_diffusion_explicit.h:71
DConvectionDiffusionExplicit()
Definition: d_convection_diffusion_explicit.h:190
void Initialize(const ProcessInfo &rCurrentProcessInfo) override
Definition: d_convection_diffusion_explicit.cpp:120
std::vector< Dof< double >::Pointer > DofsVectorType
Definition: d_convection_diffusion_explicit.h:77
std::vector< std::size_t > EquationIdVectorType
Definition: d_convection_diffusion_explicit.h:76
void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: d_convection_diffusion_explicit.cpp:68
std::size_t IndexType
Definition: d_convection_diffusion_explicit.h:75
GeometryData::IntegrationMethod IntegrationMethod
Definition: d_convection_diffusion_explicit.h:78
QSConvectionDiffusionExplicit< TDim, TNumNodes > BaseType
Definition: d_convection_diffusion_explicit.h:68
std::string Info() const override
Turn back information as a string.
Definition: d_convection_diffusion_explicit.h:150
void FinalizeSolutionStep(const ProcessInfo &rCurrentProcessInfo) override
Definition: d_convection_diffusion_explicit.cpp:135
Element::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, Properties::Pointer pProperties) const override
Definition: d_convection_diffusion_explicit.cpp:45
Matrix MatrixType
Definition: d_convection_diffusion_explicit.h:74
void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: d_convection_diffusion_explicit.cpp:82
void AddExplicitContribution(const ProcessInfo &rCurrentProcessInfo) override
Definition: d_convection_diffusion_explicit.cpp:95
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(DConvectionDiffusionExplicit)
Pointer definition of DConvectionDiffusionExplicit.
void Calculate(const Variable< double > &rVariable, double &Output, const ProcessInfo &rCurrentProcessInfo) override
Definition: d_convection_diffusion_explicit.cpp:167
Base class for all Elements.
Definition: element.h:60
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: element.h:1135
std::size_t IndexType
Definition: flags.h:74
IntegrationMethod
Definition: geometry_data.h:76
Geometry base class.
Definition: geometry.h:71
IndexType Id() const
Definition: indexed_object.h:107
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
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
This element solves the convection-diffusion equation, stabilized with algebraic subgrid scale or ort...
Definition: qs_convection_diffusion_explicit.h:71
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
#define KRATOS_SERIALIZE_SAVE_BASE_CLASS(Serializer, BaseType)
Definition: define.h:812
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
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
ProcessInfo
Definition: edgebased_PureConvection.py:116
def load(f)
Definition: ode_solve.py:307