37 #ifndef STABILIZATION_HPP 38 #define STABILIZATION_HPP 1
40 #include <lifev/core/util/Factory.hpp> 41 #include <lifev/core/util/FactorySingleton.hpp> 43 #include <lifev/core/fem/FESpace.hpp> 44 #include <lifev/core/fem/ReferenceFE.hpp> 46 #include <lifev/eta/fem/ETFESpace.hpp> 48 #include <lifev/core/array/MatrixEpetra.hpp> 49 #include <lifev/core/filter/ExporterHDF5.hpp> 51 #include <lifev/navier_stokes_blocks/solver/FastAssemblerNS.hpp>
FESpace< mesh_Type, map_Type > fespace_Type
std::shared_ptr< ETFESpace_velocity > ETFESpacePtr_velocity
std::shared_ptr< vector_Type > vectorPtr_Type
FESpace - Short description here please!
ETFESpace< mesh_Type, map_Type, 3, 1 > ETFESpace_pressure
virtual void setAlpha(const Real &alpha)=0
Set the alpha coefficient of the BDF scheme used.
void importFromHDF5(std::string const &fileName, std::string const &matrixName="matrix")
Read a matrix from a HDF5 (.h5) file.
virtual void buildGraphs()
Build the graphs of each single block.
virtual void apply_matrix(const vector_Type &, const vector_Type &, const vector_Type &)
Updates the system matrix in Navier-Stokes simulations in fixed coordinates.
std::shared_ptr< ETFESpace_pressure > ETFESpacePtr_pressure
virtual void setETvelocitySpace(const ETFESpacePtr_velocity &velocityEta_fespace)=0
Set Expression Template FE space for velocity.
virtual void setExportFineScaleVelocity(ExporterHDF5< mesh_Type > &, const int &)
Set if the user wants to export the fine scale component.
virtual matrixPtr_Type const & block_01() const =0
Get block01 of the stabilization matrix.
RegionMesh< LinearTetra > mesh_Type
FactorySingleton< Factory< Stabilization, std::string > > StabilizationFactory
virtual void setUseODEfineScale(const bool &)
Set if using dynamic fine scale model.
std::shared_ptr< matrix_Type > matrixPtr_Type
virtual void apply_matrix(const vector_Type &, const vector_Type &, const vector_Type &, const vector_Type &)
Updates the jacobian matrix in Navier-Stokes simulations in ALE coordinates.
virtual void setFastAssembler(std::shared_ptr< FastAssemblerNS > &)
Set if using the fast assembler.
Epetra_Import const & importer()
Getter for the Epetra_Import.
virtual void setViscosity(const Real &viscosity)=0
Set the fluid dynamic viscosity.
virtual void apply_vector(vectorPtr_Type &, vectorPtr_Type &, const vector_Type &, const vector_Type &, const vector_Type &, const vector_Type &)
Adds to the residual the contribution coming from the SUPG stabilization.
virtual void apply_matrix(const vector_Type &, const vector_Type &)
Updates the system matrix in Navier-Stokes simulations in ALE coordinates.
virtual std::string label()=0
Get name of stabilization being used.
virtual void setVelocitySpace(fespacePtr_Type velocityFESpace)=0
Set FE space for velocity.
virtual void updateODEfineScale(const vectorPtr_Type &, const vectorPtr_Type &)
Updates the fine scale component.
virtual void computeFineScales(const vectorPtr_Type &, const vectorPtr_Type &, const vectorPtr_Type &)
Compute the fine component of the velocity and pressure.
virtual void setupODEfineScale()
Setup of the fine scale component.
virtual void setBDForder(const Real &bdfOrder)=0
Set the bdf order.
MatrixEpetra< Real > matrix_Type
virtual matrixPtr_Type const & block_10() const =0
Get block10 of the stabilization matrix.
virtual void setPressureSpace(fespacePtr_Type pressureFESpace)=0
Set Expression Template FE space for velocity.
virtual void setUseGraph(const bool &)
Set if using the graph.
ETFESpace< mesh_Type, map_Type, 3, 3 > ETFESpace_velocity
double Real
Generic real data.
virtual matrixPtr_Type const & block_00() const =0
Get block00 of the stabilization matrix.
virtual void setDensity(const Real &density)=0
Set the fluid density.
virtual void apply_matrix(const vector_Type &)
Updates the system matrix in Navier-Stokes simulations in fixed coordinates.
virtual void setConstant(const int &value)=0
Set the constant C_I for the stabilization.
virtual void setETpressureSpace(const ETFESpacePtr_pressure &pressureEta_fespace)=0
Set Expression Template FE space for velocity.
virtual void apply_vector(vectorPtr_Type &, vectorPtr_Type &, const vector_Type &, const vector_Type &)
Adds to the right hand side the contribution coming from the SUPG stabilization.
virtual void apply_vector(vectorPtr_Type &, vectorPtr_Type &, const vector_Type &, const vector_Type &, const vector_Type &)
Adds to the residual the contribution coming from the SUPG stabilization.
virtual void setTimeStep(const Real ×tep)=0
Set the time step size.
std::shared_ptr< fespace_Type > fespacePtr_Type
virtual void computeFineScalesForVisualization(const vectorPtr_Type &, const vectorPtr_Type &, const vectorPtr_Type &)
Compute the fine component of the velocity and pressure for visualization.
class ETFESpace A light, templated version of the FESpace
virtual matrixPtr_Type const & block_11() const =0
Get block11 of the stabilization matrix.
virtual void setCommunicator(std::shared_ptr< Epetra_Comm > comm)=0
Set the Epetra communicator.
virtual void updateODEfineScale(const vectorPtr_Type &, const vectorPtr_Type &, const vectorPtr_Type &)
Updates the fine scale component.