HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
Loading...
Searching...
No Matches
Public Types | Public Member Functions | List of all members
HArDCore3D::NavierStokes Struct Reference

Class to assemble and solve a Navier-Stokes problem. More...

#include <sddr-navier-stokes.hpp>

Public Types

typedef Eigen::SparseMatrix< doubleSystemMatrixType
 
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< StokesNormscomputeStokesNorms (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, StokesNormscomputeContinuousErrorsNorms (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 SXGradsxGrad () const
 Returns the space SXGrad.
 
const SXCurlsxCurl () const
 Returns the space SXCurl.
 
const SXDivsxDiv () const
 Returns the space XDiv.
 
const SystemMatrixTypesystemMatrix () const
 Returns the linear system matrix.
 
SystemMatrixTypesystemMatrix ()
 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 SystemMatrixTypescMatrix () 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 doublestabilizationParameter () const
 Returns the stabilization parameter.
 
doublestabilizationParameter ()
 Returns the stabilization parameter.
 
const booliterativeSolver () const
 Returns the selection of iterative solver.
 
booliterativeSolver ()
 Returns the selection of iterative solver.
 
doubleReynolds ()
 Returns the reynolds number.
 

Detailed Description

Class to assemble and solve a Navier-Stokes problem.

Member Typedef Documentation

◆ ForcingTermType

typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::NavierStokes::ForcingTermType

◆ PressureGradientType

typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::NavierStokes::PressureGradientType

◆ PressureType

typedef std::function<double(const Eigen::Vector3d &)> HArDCore3D::NavierStokes::PressureType

◆ SystemMatrixType

◆ VelocityType

typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::NavierStokes::VelocityType

◆ ViscosityType

◆ VorticityType

typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::NavierStokes::VorticityType

Constructor & Destructor Documentation

◆ NavierStokes()

NavierStokes::NavierStokes ( const DDRCore ddrcore,
const BoundaryConditions BC,
bool  use_threads,
std::ostream &  output = std::cout 
)

Constructor.

Parameters
ddrcoreCore for the DDR space sequence
BCType of BC
use_threadsTrue for parallel execution, false for sequential execution
outputOutput stream to print status messages

Member Function Documentation

◆ assembleLinearSystem()

double NavierStokes::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)

Parameters
XnState in the Newton iteration
reacCoefficient of mass matrix for pseudo-temporal iterations

◆ computeContinuousErrorsNorms()

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)

Parameters
vThe vector, contains both velocity and pressure
uContinuous velocity
curl_uCurl of continuous velocity
pContinuous pressure
grad_pGrad of continuous pressure

◆ computeFlux()

double NavierStokes::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.

Parameters
uDiscrete vector (in Xcurl) of the velocity
surfaceSurface, described as a convex hull, across which the flux is computed

◆ computeLinearComponents()

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.

Parameters
fForcing term
uBoundary condition on velocity
omegaBoundary condition on vorticity
reynoldsReynolds number for BC on vorticity

◆ computeStokesNorms()

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.

Parameters
list_dofsList of vectors representing velocities and pressures

◆ iterativeSolver() [1/2]

bool & HArDCore3D::NavierStokes::iterativeSolver ( )
inline

Returns the selection of iterative solver.

◆ iterativeSolver() [2/2]

const bool & HArDCore3D::NavierStokes::iterativeSolver ( ) const
inline

Returns the selection of iterative solver.

◆ nDOFs_dir_edge_p()

size_t HArDCore3D::NavierStokes::nDOFs_dir_edge_p ( ) const
inline

Number of DOFs on Dirichlet edges for the pressure.

◆ nDOFs_dir_edge_u()

size_t HArDCore3D::NavierStokes::nDOFs_dir_edge_u ( ) const
inline

Number of DOFs on Dirichlet edges for the velocity.

◆ nDOFs_dir_face_p()

size_t HArDCore3D::NavierStokes::nDOFs_dir_face_p ( ) const
inline

Number of DOFs on Dirichlet faces for the pressure.

◆ nDOFs_dir_face_u()

size_t HArDCore3D::NavierStokes::nDOFs_dir_face_u ( ) const
inline

Number of DOFs on Dirichlet faces for the velocity.

◆ nDOFs_dir_p()

size_t HArDCore3D::NavierStokes::nDOFs_dir_p ( ) const
inline

Total number of Dirichlet DOFs for the pressure.

◆ nDOFs_dir_u()

size_t HArDCore3D::NavierStokes::nDOFs_dir_u ( ) const
inline

Total number of Dirichlet DOFs for the velocity.

◆ nDOFs_dir_vertex_p()

size_t HArDCore3D::NavierStokes::nDOFs_dir_vertex_p ( ) const
inline

Number of DOFs on Dirichlet vertices for the pressure.

◆ nDOFs_up()

size_t HArDCore3D::NavierStokes::nDOFs_up ( ) const
inline

Number of DOFs for velocity and pressure.

◆ newtonIncrement()

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)

Parameters
solverLinear solver to use
navier_scalingScaling of the Navier term
relaxRelaxation parameter
norm_dXNorm of the increment dX (before applying relaxation)

◆ nGLUKN()

size_t HArDCore3D::NavierStokes::nGLUKN ( ) const
inline

Number of globally coupled unknowns for velocity, pressure, Lagrange multiplier.

◆ nGLUKN_p()

size_t HArDCore3D::NavierStokes::nGLUKN_p ( ) const
inline

Number of globally coupled unknowns for the pressure (statically condensed unknowns removed)

◆ nGLUKN_u()

size_t HArDCore3D::NavierStokes::nGLUKN_u ( ) const
inline

Number of globally coupled unknowns for the velocity (statically condensed unknowns removed)

◆ normGlobalRHS()

double NavierStokes::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)

Parameters
locRHSList of local contributions

◆ nSCUKN()

size_t HArDCore3D::NavierStokes::nSCUKN ( ) const
inline

Number of statically condensed unknowns (both velocity and pressure)

◆ nSCUKN_p()

size_t HArDCore3D::NavierStokes::nSCUKN_p ( ) const
inline

Number of statically condensed unknowns (pressure)

◆ nSCUKN_u()

size_t HArDCore3D::NavierStokes::nSCUKN_u ( ) const
inline

Number of statically condensed unknowns (velocity)

◆ nUKN_lambda()

size_t HArDCore3D::NavierStokes::nUKN_lambda ( ) const
inline

Number of Lagrange multiplier unknowns (1 if BC=N to impose \int p=0, 0 otherwise)

◆ nUKN_p()

size_t HArDCore3D::NavierStokes::nUKN_p ( ) const
inline

Number of pressure unknowns (statically condensed and global)

◆ nUKN_u()

size_t HArDCore3D::NavierStokes::nUKN_u ( ) const
inline

Number of velocity unknowns (statically condensed and global, but no Dirichlet)

◆ nUKN_uplambda()

size_t HArDCore3D::NavierStokes::nUKN_uplambda ( ) const
inline

Number of unknowns (velocity-pressure-Lagrange multiplier, including statically condensed and globally coupled unknowns)

◆ Reynolds()

double & HArDCore3D::NavierStokes::Reynolds ( )
inline

Returns the reynolds number.

◆ rhsSourceNaturalBC()

std::vector< Eigen::VectorXd > & HArDCore3D::NavierStokes::rhsSourceNaturalBC ( )
inline

Returns local RHS for Stokes problem (useful to compute relative residual)

◆ scMatrix()

const SystemMatrixType & HArDCore3D::NavierStokes::scMatrix ( ) const
inline

Returns the static condensation recovery operator.

◆ scVector()

Eigen::VectorXd & HArDCore3D::NavierStokes::scVector ( )
inline

Returns the static condensation rhs.

◆ stabilizationParameter() [1/2]

double & HArDCore3D::NavierStokes::stabilizationParameter ( )
inline

Returns the stabilization parameter.

◆ stabilizationParameter() [2/2]

const double & HArDCore3D::NavierStokes::stabilizationParameter ( ) const
inline

Returns the stabilization parameter.

◆ sxCurl()

const SXCurl & HArDCore3D::NavierStokes::sxCurl ( ) const
inline

Returns the space SXCurl.

◆ sxDiv()

const SXDiv & HArDCore3D::NavierStokes::sxDiv ( ) const
inline

Returns the space XDiv.

◆ sxGrad()

const SXGrad & HArDCore3D::NavierStokes::sxGrad ( ) const
inline

Returns the space SXGrad.

◆ systemMatrix() [1/2]

SystemMatrixType & HArDCore3D::NavierStokes::systemMatrix ( )
inline

Returns the linear system matrix.

◆ systemMatrix() [2/2]

const SystemMatrixType & HArDCore3D::NavierStokes::systemMatrix ( ) const
inline

Returns the linear system matrix.

◆ systemVector() [1/2]

Eigen::VectorXd & HArDCore3D::NavierStokes::systemVector ( )
inline

Returns the linear system right-hand side vector.

◆ systemVector() [2/2]

const Eigen::VectorXd & HArDCore3D::NavierStokes::systemVector ( ) const
inline

Returns the linear system right-hand side vector.

◆ vectorDirBC()

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

Parameters
uBC on velocity
pBC on pressure

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