|
HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
|
Class to assemble and solve a Navier-Stokes problem. More...
#include <sddr-navier-stokes.hpp>
Public Types | |
| typedef Eigen::SparseMatrix< double > | SystemMatrixType |
| typedef std::function< Eigen::Vector3d(const Eigen::Vector3d &)> | ForcingTermType |
| typedef std::function< Eigen::Vector3d(const Eigen::Vector3d &)> | VelocityType |
| typedef std::function< Eigen::Vector3d(const Eigen::Vector3d &)> | VorticityType |
| typedef std::function< double(const Eigen::Vector3d &)> | PressureType |
| typedef std::function< Eigen::Vector3d(const Eigen::Vector3d &)> | PressureGradientType |
| typedef IntegralWeight | ViscosityType |
Public Member Functions | |
| NavierStokes (const DDRCore &ddrcore, const BoundaryConditions &BC, bool use_threads, std::ostream &output=std::cout) | |
| Constructor. | |
| Eigen::VectorXd | vectorDirBC (const VelocityType &u, const PressureType &p) |
| Create a vector of Dirichlet boundary conditions. | |
| void | computeLinearComponents (const ForcingTermType &f, const VelocityType &u, const VorticityType &omega, const double &reynolds) |
| Create the linear part of the system ( curl-curl velocity matrix, vel-grad pressure matrix, and RHS). Needs to be called before assembleLinearSystem. | |
| double | assembleLinearSystem (const Eigen::VectorXd &Xn, const double &reac) |
| Assemble the global system at each Newton iteration, and returns the residual (norm of F(Xn)-b) | |
| double | normGlobalRHS (const std::vector< Eigen::VectorXd > &locRHS) |
| Computes the norm of a vector (in Xcurl-Xgrad space) from the element-wise contributions to that vector, excluding the Dirichlet components (only the UKN are taken into account in the norm) | |
| Eigen::VectorXd | newtonIncrement (LinearSolver< SystemMatrixType > &solver, const double &navier_scaling, const double &relax, double &norm_dX) |
| Solve the linear system and computes the Newton increment (taking into account relaxation) | |
| std::vector< StokesNorms > | computeStokesNorms (const std::vector< Eigen::VectorXd > &list_dofs) const |
| Compute the discrete L2 norms, for a family of Eigen::VectorXd representing velocities & pressures, of u, curl u, p and grad p. | |
| std::pair< StokesNorms, StokesNorms > | computeContinuousErrorsNorms (const Eigen::VectorXd &v, const VelocityType &u, const VorticityType &curl_u, const PressureType &p, const PressureGradientType &grad_p) const |
| Compute the continuous Hcurl errors in u and Hgrad error in p (1st component), and same with continuous norms (second component) | |
| double | computeFlux (const Eigen::VectorXd &u, const Hull &surface) const |
| Compute the flux (based on the polynomial reconstructions) for a discrete velocity across a flat convex hull. | |
| size_t | nDOFs_up () const |
| Number of DOFs for velocity and pressure. | |
| size_t | nDOFs_dir_edge_u () const |
| Number of DOFs on Dirichlet edges for the velocity. | |
| size_t | nDOFs_dir_face_u () const |
| Number of DOFs on Dirichlet faces for the velocity. | |
| size_t | nDOFs_dir_u () const |
| Total number of Dirichlet DOFs for the velocity. | |
| size_t | nDOFs_dir_vertex_p () const |
| Number of DOFs on Dirichlet vertices for the pressure. | |
| size_t | nDOFs_dir_edge_p () const |
| Number of DOFs on Dirichlet edges for the pressure. | |
| size_t | nDOFs_dir_face_p () const |
| Number of DOFs on Dirichlet faces for the pressure. | |
| size_t | nDOFs_dir_p () const |
| Total number of Dirichlet DOFs for the pressure. | |
| size_t | nSCUKN_u () const |
| Number of statically condensed unknowns (velocity) | |
| size_t | nSCUKN_p () const |
| Number of statically condensed unknowns (pressure) | |
| size_t | nSCUKN () const |
| Number of statically condensed unknowns (both velocity and pressure) | |
| size_t | nUKN_u () const |
| Number of velocity unknowns (statically condensed and global, but no Dirichlet) | |
| size_t | nUKN_p () const |
| Number of pressure unknowns (statically condensed and global) | |
| size_t | nUKN_lambda () const |
| Number of Lagrange multiplier unknowns (1 if BC=N to impose \int p=0, 0 otherwise) | |
| size_t | nUKN_uplambda () const |
| Number of unknowns (velocity-pressure-Lagrange multiplier, including statically condensed and globally coupled unknowns) | |
| size_t | nGLUKN_u () const |
| Number of globally coupled unknowns for the velocity (statically condensed unknowns removed) | |
| size_t | nGLUKN_p () const |
| Number of globally coupled unknowns for the pressure (statically condensed unknowns removed) | |
| size_t | nGLUKN () const |
| Number of globally coupled unknowns for velocity, pressure, Lagrange multiplier. | |
| const SXGrad & | sxGrad () const |
| Returns the space SXGrad. | |
| const SXCurl & | sxCurl () const |
| Returns the space SXCurl. | |
| const SXDiv & | sxDiv () const |
| Returns the space XDiv. | |
| const SystemMatrixType & | systemMatrix () const |
| Returns the linear system matrix. | |
| SystemMatrixType & | systemMatrix () |
| 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 SystemMatrixType & | scMatrix () const |
| Returns the static condensation recovery operator. | |
| Eigen::VectorXd & | scVector () |
| Returns the static condensation rhs. | |
| std::vector< Eigen::VectorXd > & | rhsSourceNaturalBC () |
| Returns local RHS for Stokes problem (useful to compute relative residual) | |
| const double & | stabilizationParameter () const |
| Returns the stabilization parameter. | |
| double & | stabilizationParameter () |
| Returns the stabilization parameter. | |
| const bool & | iterativeSolver () const |
| Returns the selection of iterative solver. | |
| bool & | iterativeSolver () |
| Returns the selection of iterative solver. | |
| double & | Reynolds () |
| Returns the reynolds number. | |
Class to assemble and solve a Navier-Stokes problem.
| typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::NavierStokes::ForcingTermType |
| typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::NavierStokes::PressureGradientType |
| typedef std::function<double(const Eigen::Vector3d &)> HArDCore3D::NavierStokes::PressureType |
| typedef Eigen::SparseMatrix<double> HArDCore3D::NavierStokes::SystemMatrixType |
| typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::NavierStokes::VelocityType |
| typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::NavierStokes::VorticityType |
| NavierStokes::NavierStokes | ( | const DDRCore & | ddrcore, |
| const BoundaryConditions & | BC, | ||
| bool | use_threads, | ||
| std::ostream & | output = std::cout |
||
| ) |
Constructor.
| ddrcore | Core for the DDR space sequence |
| BC | Type of BC |
| use_threads | True for parallel execution, false for sequential execution |
| output | Output stream to print status messages |
Assemble the global system at each Newton iteration, and returns the residual (norm of F(Xn)-b)
| Xn | State in the Newton iteration |
| reac | Coefficient of mass matrix for pseudo-temporal iterations |
| std::pair< StokesNorms, StokesNorms > NavierStokes::computeContinuousErrorsNorms | ( | const Eigen::VectorXd & | v, |
| const VelocityType & | u, | ||
| const VorticityType & | curl_u, | ||
| const PressureType & | p, | ||
| const PressureGradientType & | grad_p | ||
| ) | const |
Compute the continuous Hcurl errors in u and Hgrad error in p (1st component), and same with continuous norms (second component)
| v | The vector, contains both velocity and pressure |
| u | Continuous velocity |
| curl_u | Curl of continuous velocity |
| p | Continuous pressure |
| grad_p | Grad of continuous pressure |
Compute the flux (based on the polynomial reconstructions) for a discrete velocity across a flat convex hull.
| u | Discrete vector (in Xcurl) of the velocity |
| surface | Surface, described as a convex hull, across which the flux is computed |
| void NavierStokes::computeLinearComponents | ( | const ForcingTermType & | f, |
| const VelocityType & | u, | ||
| const VorticityType & | omega, | ||
| const double & | reynolds | ||
| ) |
Create the linear part of the system ( curl-curl velocity matrix, vel-grad pressure matrix, and RHS). Needs to be called before assembleLinearSystem.
| f | Forcing term |
| u | Boundary condition on velocity |
| omega | Boundary condition on vorticity |
| reynolds | Reynolds number for BC on vorticity |
| std::vector< StokesNorms > NavierStokes::computeStokesNorms | ( | const std::vector< Eigen::VectorXd > & | list_dofs | ) | const |
Compute the discrete L2 norms, for a family of Eigen::VectorXd representing velocities & pressures, of u, curl u, p and grad p.
| list_dofs | List of vectors representing velocities and pressures |
|
inline |
Returns the selection of iterative solver.
Returns the selection of iterative solver.
|
inline |
Number of DOFs on Dirichlet edges for the pressure.
|
inline |
Number of DOFs on Dirichlet edges for the velocity.
|
inline |
Number of DOFs on Dirichlet faces for the pressure.
|
inline |
Number of DOFs on Dirichlet faces for the velocity.
|
inline |
Total number of Dirichlet DOFs for the pressure.
|
inline |
Total number of Dirichlet DOFs for the velocity.
|
inline |
Number of DOFs on Dirichlet vertices for the pressure.
|
inline |
Number of DOFs for velocity and pressure.
| Eigen::VectorXd NavierStokes::newtonIncrement | ( | LinearSolver< SystemMatrixType > & | solver, |
| const double & | navier_scaling, | ||
| const double & | relax, | ||
| double & | norm_dX | ||
| ) |
Solve the linear system and computes the Newton increment (taking into account relaxation)
| solver | Linear solver to use |
| navier_scaling | Scaling of the Navier term |
| relax | Relaxation parameter |
| norm_dX | Norm of the increment dX (before applying relaxation) |
|
inline |
Number of globally coupled unknowns for velocity, pressure, Lagrange multiplier.
|
inline |
Number of globally coupled unknowns for the pressure (statically condensed unknowns removed)
|
inline |
Number of globally coupled unknowns for the velocity (statically condensed unknowns removed)
Computes the norm of a vector (in Xcurl-Xgrad space) from the element-wise contributions to that vector, excluding the Dirichlet components (only the UKN are taken into account in the norm)
| locRHS | List of local contributions |
|
inline |
Number of statically condensed unknowns (both velocity and pressure)
|
inline |
Number of statically condensed unknowns (pressure)
|
inline |
Number of statically condensed unknowns (velocity)
|
inline |
Number of Lagrange multiplier unknowns (1 if BC=N to impose \int p=0, 0 otherwise)
|
inline |
Number of pressure unknowns (statically condensed and global)
|
inline |
Number of velocity unknowns (statically condensed and global, but no Dirichlet)
|
inline |
Number of unknowns (velocity-pressure-Lagrange multiplier, including statically condensed and globally coupled unknowns)
|
inline |
Returns the reynolds number.
|
inline |
Returns local RHS for Stokes problem (useful to compute relative residual)
|
inline |
Returns the static condensation recovery operator.
|
inline |
Returns the static condensation rhs.
|
inline |
Returns the stabilization parameter.
Returns the stabilization parameter.
|
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.
| Eigen::VectorXd NavierStokes::vectorDirBC | ( | const VelocityType & | u, |
| const PressureType & | p | ||
| ) |
Create a vector of Dirichlet boundary conditions.
The vector has dimension velocity-pressure and zero everywhere except for Dirichlet DOFs on the velocity and pressure
| u | BC on velocity |
| p | BC on pressure |