|
HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
|
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 VHHOSpace & | vhhospace () const |
| Returns the velocity space. | |
| const GlobalDOFSpace & | pspace () const |
| Returns the pressure space. | |
| const GlobalDOFSpace & | phispace () const |
| Returns the pressure space. | |
| const Mesh & | mesh () const |
| Returns the mesh. | |
| const systemMatrixType & | systemMatrixHHO () const |
| Returns the linear system matrix. | |
| systemMatrixType & | systemMatrixHHO () |
| 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 systemMatrixType & | scMatrix () const |
| Returns the static condensation recovery operator. | |
| Eigen::VectorXd & | scVector () |
| Returns the static condensation rhs. | |
| const systemMatrixType & | systemMatrixDG () const |
| Returns the linear system matrix. | |
| systemMatrixType & | systemMatrixDG () |
| 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. | |
| std::vector< double > | computePhysicalEnergies (const Eigen::VectorXd &uprho, const MomentumForcingTermType &f, const ViscosityType &mu) const |
| UTILS. | |
| 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 °ree, const int &time_it, const double &reynolds, const std::vector< std::string > &plot_variables, const std::string &output_path) |
| Write solution to vut file. | |
| void | writeTimeStepData (const std::string &filename, const int &time_it, const double &time, const double &dt, const double &L2_u_err, const double &H1_u_err, const double &Energy_err, const double &L2_rho_err, const double &cf_rho_err, const double &upw_rho_err, const double &E_kinetic, const double &E_potential, const double &E_dissipation) const |
| Writing data to txt file. | |
| void | writeFileHeader (const std::string &filename, const std::string &solution_name, const int &solution_number, const std::string &mesh_name, const size_t °ree, const std::unique_ptr< Mesh > &mesh_ptr, const double &viscosity) const |
| Writing file header with simulation metadata. | |
| 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 | |
Class for Navier-Stokes model.
The global unknowns are in the following order:
| typedef std::function<double(const VectorRd &)> HArDCore2D::VarDensityNavierStokes::CompressibilityForcingTermType |
| typedef std::function<double(const VectorRd &)> HArDCore2D::VarDensityNavierStokes::DensityForcingTermType |
| typedef std::function<VectorRd(const VectorRd &)> HArDCore2D::VarDensityNavierStokes::MomentumForcingTermType |
| typedef std::function<VectorRd(const VectorRd &)> HArDCore2D::VarDensityNavierStokes::PressureGradientType |
| typedef std::function<double(const VectorRd &)> HArDCore2D::VarDensityNavierStokes::PressureType |
| typedef Eigen::SparseMatrix<double> HArDCore2D::VarDensityNavierStokes::systemMatrixType |
| typedef std::function<MatrixRd(const VectorRd &)> HArDCore2D::VarDensityNavierStokes::VelocityGradientType |
| typedef std::function<VectorRd(const VectorRd &)> HArDCore2D::VarDensityNavierStokes::VelocityType |
| typedef std::function<double(const VectorRd &)> HArDCore2D::VarDensityNavierStokes::VolumicFractionType |
| 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.
| vhho_space | Velocity space and operators |
| p_space | Space for the pressure |
| phi_space | Space for the volumic fraction |
| BC | Boundary conditions |
| use_threads | True for parallel execution, false for sequential execution |
| stokes | True if we want to assemble the Stokes system without convective terms |
| output | Output stream to print status messages |
| void VarDensityNavierStokes::assembleLinearDGSystem | ( | const CompressibilityForcingTermType & | h, |
| const Eigen::VectorXd & | uprho, | ||
| const VolumicFractionType & | rho, | ||
| const double & | time_step | ||
| ) |
DG ADVECTION SYSTEM.
| 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.
| f | Forcing term |
| g | Forcing term for the divergence equation |
| rhog | denisty on the boundary |
| mu | Viscosity |
| u0 | solution at previous time step |
| uold | solution at previous Newton iteration |
| time_step | time step |
| 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.
| 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.
| list_dofs | The list of vectors representing the dofs |
| std::vector< double > VarDensityNavierStokes::computePhysicalEnergies | ( | const Eigen::VectorXd & | uprho, |
| const MomentumForcingTermType & | f, | ||
| const ViscosityType & | mu | ||
| ) | const |
UTILS.
| uprho | Solution vector |
| f | Momentum forcing term |
| mu | Viscosity |
| std::vector< double > VarDensityNavierStokes::computePressureL2Norm | ( | const std::vector< Eigen::VectorXd > & | list_dofs | ) | const |
Compute the discrete L2 norm of the pressure.
| list_dofs | The list of vectors representing the dofs |
|
inline |
Returns the dimension of pressure space.
|
inline |
Returns the dimension of velocity space.
|
inline |
Returns the dimension of pressure space.
| double HArDCore2D::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
| Eigen::VectorXd VarDensityNavierStokes::imposingBC | ( | const VelocityType & | u, |
| const Eigen::VectorXd & | uph | ||
| ) | const |
impose boundary conditions
| 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.
| u | Continuous velocity to interpolate |
| p | Continuous pressure to interpolate |
| doe_cell | The optional degre of cell quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used. |
| doe_face | The optional degre of face quadrature rules to compute the interpolate. If negative, then 2*degree()+3 will be used. |
|
inline |
Returns the Stokes status.
|
inline |
Returns the mesh.
| 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
|
inline |
Returns the local number of pressure statically condensed DOFs.
|
inline |
Returns the local number of velocity statically condensed DOFs.
|
inline |
Returns the number of Dirichlet DOFs.
|
inline |
Returns the number of DOFs after SC and with Lagrange multiplier, but before eliminating Dirichlet DOFs.
|
inline |
Returns the number of statically condensed DOFs.
|
inline |
Returns the number of pressure statically condensed DOFs.
|
inline |
Returns the number of velocity statically condensed DOFs.
|
inline |
Returns the pressure space.
| Eigen::VectorXd VarDensityNavierStokes::pressureVertexValues | ( | const Eigen::VectorXd & | p | ) | const |
Create vertex values for the pressure (from the element values), for plotting.
| p | Vector of pressure DOFs |
|
inline |
Returns the pressure space.
|
inline |
Returns the static condensation recovery operator.
|
inline |
Returns the static condensation rhs.
|
inline |
|
inline |
Returns the size of the final system with Lagrange multiplier, after application of SC and removal of Dirichlet BCs.
|
inline |
Returns the stabilization parameter.
|
inline |
Returns the stabilization parameter (scaling)
|
inline |
Returns the linear system matrix.
|
inline |
Returns the linear system matrix.
|
inline |
Returns the linear system matrix.
|
inline |
Returns the linear system matrix.
|
inline |
Returns the linear system right-hand side vector.
|
inline |
Returns the linear system right-hand side vector.
|
inline |
Returns the linear system right-hand side vector.
|
inline |
Returns the linear system right-hand side vector.
|
inline |
Returns the velocity space.
| Eigen::VectorXd VarDensityNavierStokes::volumicFractionVertexValues | ( | const Eigen::VectorXd & | phi | ) | const |
Create vertex values for volumic fraction.
Creating vertex values for volumic fraction.
| phi | Vector of volumic fraction DOFs |
| void VarDensityNavierStokes::writeFileHeader | ( | const std::string & | filename, |
| const std::string & | solution_name, | ||
| const int & | solution_number, | ||
| const std::string & | mesh_name, | ||
| const size_t & | degree, | ||
| const std::unique_ptr< Mesh > & | mesh_ptr, | ||
| const double & | viscosity | ||
| ) | const |
Writing file header with simulation metadata.
| void VarDensityNavierStokes::writeTimeStepData | ( | const std::string & | filename, |
| const int & | time_it, | ||
| const double & | time, | ||
| const double & | dt, | ||
| const double & | L2_u_err, | ||
| const double & | H1_u_err, | ||
| const double & | Energy_err, | ||
| const double & | L2_rho_err, | ||
| const double & | cf_rho_err, | ||
| const double & | upw_rho_err, | ||
| const double & | E_kinetic, | ||
| const double & | E_potential, | ||
| const double & | E_dissipation | ||
| ) | const |
Writing data to txt file.
| 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, | ||
| const double & | reynolds, | ||
| const std::vector< std::string > & | plot_variables, | ||
| const std::string & | output_path | ||
| ) |
Write solution to vut file.
Writing solution to vtu file.