HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
Loading...
Searching...
No Matches
Public Types | Public Member Functions | List of all members
HArDCore2D::VarDensityNavierStokes Struct Reference

Class for Navier-Stokes model. More...

#include <hho-nsa.hpp>

Public Types

typedef Eigen::SparseMatrix< double > systemMatrixType
 
typedef std::function< VectorRd(const VectorRd &)> MomentumForcingTermType
 
typedef std::function< double(const VectorRd &)> CompressibilityForcingTermType
 
typedef std::function< double(const VectorRd &)> DensityForcingTermType
 
typedef std::function< VectorRd(const VectorRd &)> VelocityType
 
typedef std::function< MatrixRd(const VectorRd &)> VelocityGradientType
 
typedef std::function< double(const VectorRd &)> PressureType
 
typedef std::function< VectorRd(const VectorRd &)> PressureGradientType
 
typedef std::function< double(const VectorRd &)> VolumicFractionType
 
typedef IntegralWeight ViscosityType
 

Public Member Functions

 VarDensityNavierStokes (const VHHOSpace &vhho_space, const GlobalDOFSpace &p_space, const GlobalDOFSpace &phi_space, const BoundaryConditions &BC, bool use_threads, bool stokes=false, std::ostream &output=std::cout)
 Constructor.
 
void assembleLinearHHOSystem (const MomentumForcingTermType &f, const CompressibilityForcingTermType &g, const VolumicFractionType &rhog, const ViscosityType &mu, const Eigen::VectorXd &u0, const Eigen::VectorXd &uold, const double &time_step)
 Assemble the global system

 
void assembleLinearDGSystem (const CompressibilityForcingTermType &h, const Eigen::VectorXd &uprho, const VolumicFractionType &rho, const double &time_step)
 DG ADVECTION SYSTEM.
 
size_t nloc_sc_u () const
 Returns the local number of velocity statically condensed DOFs.
 
size_t nloc_sc_p () const
 Returns the local number of pressure statically condensed DOFs.
 
size_t numSCDOFs_u () const
 Returns the number of velocity statically condensed DOFs.
 
size_t numSCDOFs_p () const
 Returns the number of pressure statically condensed DOFs.
 
size_t numSCDOFs () const
 Returns the number of statically condensed DOFs.
 
size_t numDirDOFs () const
 Returns the number of Dirichlet DOFs.
 
size_t dimVelocity () const
 Returns the dimension of velocity space.
 
size_t dimPressure () const
 Returns the dimension of pressure space.
 
size_t dimVolumicFraction () const
 Returns the dimension of pressure space.
 
size_t numNonSCDOFs () const
 Returns the number of DOFs after SC and with Lagrange multiplier, but before eliminating Dirichlet DOFs.
 
size_t sizeSystemHHO () const
 Returns the size of the final system with Lagrange multiplier, after application of SC and removal of Dirichlet BCs.
 
size_t sizeSystemDG () const
 
const VHHOSpacevhhospace () const
 Returns the velocity space.
 
const GlobalDOFSpacepspace () const
 Returns the pressure space.
 
const GlobalDOFSpacephispace () const
 Returns the pressure space.
 
const Meshmesh () const
 Returns the mesh.
 
const systemMatrixTypesystemMatrixHHO () const
 Returns the linear system matrix.
 
systemMatrixTypesystemMatrixHHO ()
 Returns the linear system matrix.
 
const Eigen::VectorXd & systemVector () const
 Returns the linear system right-hand side vector.
 
Eigen::VectorXd & systemVector ()
 Returns the linear system right-hand side vector.
 
const double & stabilizationParameter () const
 Returns the stabilization parameter (scaling)
 
double & stabilizationParameter ()
 Returns the stabilization parameter.
 
const systemMatrixTypescMatrix () const
 Returns the static condensation recovery operator.
 
Eigen::VectorXd & scVector ()
 Returns the static condensation rhs.
 
const systemMatrixTypesystemMatrixDG () const
 Returns the linear system matrix.
 
systemMatrixTypesystemMatrixDG ()
 Returns the linear system matrix.
 
const Eigen::VectorXd & systemVectorDG () const
 Returns the linear system right-hand side vector.
 
Eigen::VectorXd & systemVectorDG ()
 Returns the linear system right-hand side vector.
 
bool & isStokes ()
 Returns the Stokes status.
 
Eigen::VectorXd interpolate (const VelocityType &u, const PressureType &p, const VolumicFractionType &phi, const int doe_cell=-1, const int doe_face=-1) const
 Interpolates velocity and pressure.
 
std::vector< double > computeEnergyNorms (const std::vector< Eigen::VectorXd > &list_dofs) const
 Compute the discrete energy norm of a family of vectors representing the dofs.
 
std::vector< double > computePressureL2Norm (const std::vector< Eigen::VectorXd > &list_dofs) const
 Compute the discrete L2 norm of the pressure.
 
std::vector< std::tuple< double, double, double > > computeDGNorms (const std::vector< Eigen::VectorXd > &list_dofs, const VelocityType &u) const
 Compute the discrete energy norm of a family of vectors representing the dofs.
 
Eigen::VectorXd pressureVertexValues (const Eigen::VectorXd &p) const
 Create vertex values for the pressure (from the element values), for plotting.
 
Eigen::VectorXd volumicFractionVertexValues (const Eigen::VectorXd &phi) const
 Create vertex values for volumic fraction.
 
void writeToVtuFile (const Eigen::VectorXd &upphi, const Eigen::VectorXd &upphiI, const std::unique_ptr< Mesh > &mesh_ptr, const std::string &solution_name, const std::string &mesh_name, const size_t &degree, const int &time_it)
 Write solution to vut file.
 
Eigen::VectorXd imposingBC (const VelocityType &u, const Eigen::VectorXd &uph) const
 impose boundary conditions
 
double getUpwindDensity (const size_t &iE, const size_t &iT, const Eigen::VectorXd &uprho) const
 get upwind value of the density at the face E
 
Eigen::VectorXd newtonRaphson (LinearSolver< VarDensityNavierStokes::systemMatrixType > &solver, const std::unique_ptr< Mesh > &mesh_ptr, const size_t &newton_maxit, const double &newton_tol, const Eigen::VectorXd &u_newt, const Eigen::VectorXd &un, const MomentumForcingTermType &f, const CompressibilityForcingTermType &g, const VolumicFractionType &rhog, const ViscosityType &mu, const double &time_step)
 newton method
 

Detailed Description

Class for Navier-Stokes model.

The global unknowns are in the following order:

Member Typedef Documentation

◆ CompressibilityForcingTermType

◆ DensityForcingTermType

◆ MomentumForcingTermType

◆ PressureGradientType

◆ PressureType

typedef std::function<double(const VectorRd &)> HArDCore2D::VarDensityNavierStokes::PressureType

◆ systemMatrixType

typedef Eigen::SparseMatrix<double> HArDCore2D::VarDensityNavierStokes::systemMatrixType

◆ VelocityGradientType

◆ VelocityType

◆ ViscosityType

◆ VolumicFractionType

typedef std::function<double(const VectorRd &)> HArDCore2D::VarDensityNavierStokes::VolumicFractionType

Constructor & Destructor Documentation

◆ VarDensityNavierStokes()

VarDensityNavierStokes::VarDensityNavierStokes ( const VHHOSpace vhho_space,
const GlobalDOFSpace p_space,
const GlobalDOFSpace phi_space,
const BoundaryConditions BC,
bool  use_threads,
bool  stokes = false,
std::ostream &  output = std::cout 
)

Constructor.

Parameters
vhho_spaceVelocity space and operators
p_spaceSpace for the pressure
phi_spaceSpace for the volumic fraction
BCBoundary conditions
use_threadsTrue for parallel execution, false for sequential execution
stokesTrue if we want to assemble the Stokes system without convective terms
outputOutput stream to print status messages

Member Function Documentation

◆ assembleLinearDGSystem()

void VarDensityNavierStokes::assembleLinearDGSystem ( const CompressibilityForcingTermType h,
const Eigen::VectorXd &  uprho,
const VolumicFractionType rho,
const double &  time_step 
)

DG ADVECTION SYSTEM.

◆ assembleLinearHHOSystem()

void VarDensityNavierStokes::assembleLinearHHOSystem ( const MomentumForcingTermType f,
const CompressibilityForcingTermType g,
const VolumicFractionType rhog,
const ViscosityType mu,
const Eigen::VectorXd &  u0,
const Eigen::VectorXd &  uold,
const double &  time_step 
)

Assemble the global system

HHO NAVIER-STOKES SYSTEM.

Parameters
fForcing term
gForcing term for the divergence equation
rhogdenisty on the boundary
muViscosity
u0solution at previous time step
uoldsolution at previous Newton iteration
time_steptime step

◆ computeDGNorms()

std::vector< std::tuple< double, double, double > > VarDensityNavierStokes::computeDGNorms ( const std::vector< Eigen::VectorXd > &  list_dofs,
const VelocityType u 
) const

Compute the discrete energy norm of a family of vectors representing the dofs.

◆ computeEnergyNorms()

std::vector< double > VarDensityNavierStokes::computeEnergyNorms ( const std::vector< Eigen::VectorXd > &  list_dofs) const

Compute the discrete energy norm of a family of vectors representing the dofs.

Parameters
list_dofsThe list of vectors representing the dofs

◆ computePressureL2Norm()

std::vector< double > VarDensityNavierStokes::computePressureL2Norm ( const std::vector< Eigen::VectorXd > &  list_dofs) const

Compute the discrete L2 norm of the pressure.

Parameters
list_dofsThe list of vectors representing the dofs

◆ dimPressure()

size_t HArDCore2D::VarDensityNavierStokes::dimPressure ( ) const
inline

Returns the dimension of pressure space.

◆ dimVelocity()

size_t HArDCore2D::VarDensityNavierStokes::dimVelocity ( ) const
inline

Returns the dimension of velocity space.

◆ dimVolumicFraction()

size_t HArDCore2D::VarDensityNavierStokes::dimVolumicFraction ( ) const
inline

Returns the dimension of pressure space.

◆ getUpwindDensity()

double VarDensityNavierStokes::getUpwindDensity ( const size_t &  iE,
const size_t &  iT,
const Eigen::VectorXd &  uprho 
) const

get upwind value of the density at the face E

◆ imposingBC()

Eigen::VectorXd VarDensityNavierStokes::imposingBC ( const VelocityType u,
const Eigen::VectorXd &  uph 
) const

impose boundary conditions

◆ interpolate()

Eigen::VectorXd VarDensityNavierStokes::interpolate ( const VelocityType u,
const PressureType p,
const VolumicFractionType phi,
const int  doe_cell = -1,
const int  doe_face = -1 
) const

Interpolates velocity and pressure.

INTERPOLATION.

Parameters
uContinuous velocity to interpolate
pContinuous pressure to interpolate
doe_cellThe optional degre of cell quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used.
doe_faceThe optional degre of face quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used.

◆ isStokes()

bool & HArDCore2D::VarDensityNavierStokes::isStokes ( )
inline

Returns the Stokes status.

◆ mesh()

const Mesh & HArDCore2D::VarDensityNavierStokes::mesh ( ) const
inline

Returns the mesh.

◆ newtonRaphson()

Eigen::VectorXd VarDensityNavierStokes::newtonRaphson ( LinearSolver< VarDensityNavierStokes::systemMatrixType > &  solver,
const std::unique_ptr< Mesh > &  mesh_ptr,
const size_t &  newton_maxit,
const double &  newton_tol,
const Eigen::VectorXd &  u_newt,
const Eigen::VectorXd &  un,
const MomentumForcingTermType f,
const CompressibilityForcingTermType g,
const VolumicFractionType rhog,
const ViscosityType mu,
const double &  time_step 
)

newton method

initializing parameters

◆ nloc_sc_p()

size_t HArDCore2D::VarDensityNavierStokes::nloc_sc_p ( ) const
inline

Returns the local number of pressure statically condensed DOFs.

◆ nloc_sc_u()

size_t HArDCore2D::VarDensityNavierStokes::nloc_sc_u ( ) const
inline

Returns the local number of velocity statically condensed DOFs.

◆ numDirDOFs()

size_t HArDCore2D::VarDensityNavierStokes::numDirDOFs ( ) const
inline

Returns the number of Dirichlet DOFs.

◆ numNonSCDOFs()

size_t HArDCore2D::VarDensityNavierStokes::numNonSCDOFs ( ) const
inline

Returns the number of DOFs after SC and with Lagrange multiplier, but before eliminating Dirichlet DOFs.

◆ numSCDOFs()

size_t HArDCore2D::VarDensityNavierStokes::numSCDOFs ( ) const
inline

Returns the number of statically condensed DOFs.

◆ numSCDOFs_p()

size_t HArDCore2D::VarDensityNavierStokes::numSCDOFs_p ( ) const
inline

Returns the number of pressure statically condensed DOFs.

◆ numSCDOFs_u()

size_t HArDCore2D::VarDensityNavierStokes::numSCDOFs_u ( ) const
inline

Returns the number of velocity statically condensed DOFs.

◆ phispace()

const GlobalDOFSpace & HArDCore2D::VarDensityNavierStokes::phispace ( ) const
inline

Returns the pressure space.

◆ pressureVertexValues()

Eigen::VectorXd VarDensityNavierStokes::pressureVertexValues ( const Eigen::VectorXd &  p) const

Create vertex values for the pressure (from the element values), for plotting.

Parameters
pVector of pressure DOFs

◆ pspace()

const GlobalDOFSpace & HArDCore2D::VarDensityNavierStokes::pspace ( ) const
inline

Returns the pressure space.

◆ scMatrix()

const systemMatrixType & HArDCore2D::VarDensityNavierStokes::scMatrix ( ) const
inline

Returns the static condensation recovery operator.

◆ scVector()

Eigen::VectorXd & HArDCore2D::VarDensityNavierStokes::scVector ( )
inline

Returns the static condensation rhs.

◆ sizeSystemDG()

size_t HArDCore2D::VarDensityNavierStokes::sizeSystemDG ( ) const
inline

◆ sizeSystemHHO()

size_t HArDCore2D::VarDensityNavierStokes::sizeSystemHHO ( ) const
inline

Returns the size of the final system with Lagrange multiplier, after application of SC and removal of Dirichlet BCs.

◆ stabilizationParameter() [1/2]

double & HArDCore2D::VarDensityNavierStokes::stabilizationParameter ( )
inline

Returns the stabilization parameter.

◆ stabilizationParameter() [2/2]

const double & HArDCore2D::VarDensityNavierStokes::stabilizationParameter ( ) const
inline

Returns the stabilization parameter (scaling)

◆ systemMatrixDG() [1/2]

systemMatrixType & HArDCore2D::VarDensityNavierStokes::systemMatrixDG ( )
inline

Returns the linear system matrix.

◆ systemMatrixDG() [2/2]

const systemMatrixType & HArDCore2D::VarDensityNavierStokes::systemMatrixDG ( ) const
inline

Returns the linear system matrix.

◆ systemMatrixHHO() [1/2]

systemMatrixType & HArDCore2D::VarDensityNavierStokes::systemMatrixHHO ( )
inline

Returns the linear system matrix.

◆ systemMatrixHHO() [2/2]

const systemMatrixType & HArDCore2D::VarDensityNavierStokes::systemMatrixHHO ( ) const
inline

Returns the linear system matrix.

◆ systemVector() [1/2]

Eigen::VectorXd & HArDCore2D::VarDensityNavierStokes::systemVector ( )
inline

Returns the linear system right-hand side vector.

◆ systemVector() [2/2]

const Eigen::VectorXd & HArDCore2D::VarDensityNavierStokes::systemVector ( ) const
inline

Returns the linear system right-hand side vector.

◆ systemVectorDG() [1/2]

Eigen::VectorXd & HArDCore2D::VarDensityNavierStokes::systemVectorDG ( )
inline

Returns the linear system right-hand side vector.

◆ systemVectorDG() [2/2]

const Eigen::VectorXd & HArDCore2D::VarDensityNavierStokes::systemVectorDG ( ) const
inline

Returns the linear system right-hand side vector.

◆ vhhospace()

const VHHOSpace & HArDCore2D::VarDensityNavierStokes::vhhospace ( ) const
inline

Returns the velocity space.

◆ volumicFractionVertexValues()

Eigen::VectorXd VarDensityNavierStokes::volumicFractionVertexValues ( const Eigen::VectorXd &  phi) const

Create vertex values for volumic fraction.

Creating vertex values for volumic fraction.

Parameters
phiVector of volumic fraction DOFs

◆ writeToVtuFile()

void VarDensityNavierStokes::writeToVtuFile ( const Eigen::VectorXd &  upphi,
const Eigen::VectorXd &  upphiI,
const std::unique_ptr< Mesh > &  mesh_ptr,
const std::string &  solution_name,
const std::string &  mesh_name,
const size_t &  degree,
const int &  time_it 
)

Write solution to vut file.

Writing solution to vtu file.


The documentation for this struct was generated from the following files: