Nektar++
|
Upwind scheme Riemann solver. More...
#include <UpwindSolver.h>
Static Public Member Functions | |
static SOLVER_UTILS_EXPORT RiemannSolverSharedPtr | create () |
Static Public Attributes | |
static std::string | solverName |
Protected Member Functions | |
UpwindSolver () | |
Default constructor. | |
virtual void | v_Solve (const int nDim, const Array< OneD, const Array< OneD, NekDouble > > &Fwd, const Array< OneD, const Array< OneD, NekDouble > > &Bwd, Array< OneD, Array< OneD, NekDouble > > &flux) |
Implementation of the upwind solver. | |
Protected Member Functions inherited from Nektar::SolverUtils::RiemannSolver | |
SOLVER_UTILS_EXPORT | RiemannSolver () |
void | GenerateRotationMatrices (const Array< OneD, const Array< OneD, NekDouble > > &normals) |
Generate rotation matrices for 3D expansions. | |
void | FromToRotation (Array< OneD, const NekDouble > &from, Array< OneD, const NekDouble > &to, NekDouble *mat) |
A function for creating a rotation matrix that rotates a vector from into another vector to. | |
SOLVER_UTILS_EXPORT void | rotateToNormal (const Array< OneD, const Array< OneD, NekDouble > > &inarray, const Array< OneD, const Array< OneD, NekDouble > > &normals, const Array< OneD, const Array< OneD, NekDouble > > &vecLocs, Array< OneD, Array< OneD, NekDouble > > &outarray) |
Rotate a vector field to trace normal. | |
SOLVER_UTILS_EXPORT void | rotateFromNormal (const Array< OneD, const Array< OneD, NekDouble > > &inarray, const Array< OneD, const Array< OneD, NekDouble > > &normals, const Array< OneD, const Array< OneD, NekDouble > > &vecLocs, Array< OneD, Array< OneD, NekDouble > > &outarray) |
Rotate a vector field from trace normal. | |
SOLVER_UTILS_EXPORT bool | CheckScalars (std::string name) |
Determine whether a scalar has been defined in m_scalars. | |
SOLVER_UTILS_EXPORT bool | CheckVectors (std::string name) |
Determine whether a vector has been defined in m_vectors. | |
SOLVER_UTILS_EXPORT bool | CheckParams (std::string name) |
Determine whether a parameter has been defined in m_params. | |
SOLVER_UTILS_EXPORT bool | CheckAuxScal (std::string name) |
Determine whether a scalar has been defined in m_auxScal. | |
SOLVER_UTILS_EXPORT bool | CheckAuxVec (std::string name) |
Determine whether a vector has been defined in m_auxVec. |
Additional Inherited Members | |
Public Member Functions inherited from Nektar::SolverUtils::RiemannSolver | |
SOLVER_UTILS_EXPORT void | Solve (const int nDim, const Array< OneD, const Array< OneD, NekDouble > > &Fwd, const Array< OneD, const Array< OneD, NekDouble > > &Bwd, Array< OneD, Array< OneD, NekDouble > > &flux) |
Perform the Riemann solve given the forwards and backwards spaces. | |
template<typename FuncPointerT , typename ObjectPointerT > | |
void | SetScalar (std::string name, FuncPointerT func, ObjectPointerT obj) |
void | SetScalar (std::string name, RSScalarFuncType fp) |
template<typename FuncPointerT , typename ObjectPointerT > | |
void | SetVector (std::string name, FuncPointerT func, ObjectPointerT obj) |
void | SetVector (std::string name, RSVecFuncType fp) |
template<typename FuncPointerT , typename ObjectPointerT > | |
void | SetParam (std::string name, FuncPointerT func, ObjectPointerT obj) |
void | SetParam (std::string name, RSParamFuncType fp) |
template<typename FuncPointerT , typename ObjectPointerT > | |
void | SetAuxScal (std::string name, FuncPointerT func, ObjectPointerT obj) |
template<typename FuncPointerT , typename ObjectPointerT > | |
void | SetAuxVec (std::string name, FuncPointerT func, ObjectPointerT obj) |
std::map< std::string, RSScalarFuncType > & | GetScalars () |
std::map< std::string, RSVecFuncType > & | GetVectors () |
std::map< std::string, RSParamFuncType > & | GetParams () |
Public Attributes inherited from Nektar::SolverUtils::RiemannSolver | |
int | m_spacedim |
Protected Attributes inherited from Nektar::SolverUtils::RiemannSolver | |
bool | m_requiresRotation |
Indicates whether the Riemann solver requires a rotation to be applied to the velocity fields. | |
std::map< std::string, RSScalarFuncType > | m_scalars |
Map of scalar function types. | |
std::map< std::string, RSVecFuncType > | m_vectors |
Map of vector function types. | |
std::map< std::string, RSParamFuncType > | m_params |
Map of parameter function types. | |
std::map< std::string, RSScalarFuncType > | m_auxScal |
Map of auxiliary scalar function types. | |
std::map< std::string, RSVecFuncType > | m_auxVec |
Map of auxiliary vector function types. | |
Array< OneD, Array< OneD, NekDouble > > | m_rotMat |
Rotation matrices for each trace quadrature point. | |
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > | m_rotStorage |
Rotation storage. |
Upwind scheme Riemann solver.
The upwind solver determines the flux based upon an advection velocity and trace normal . In particular, the flux for each component of the velocity field is deterined by:
Here the superscript + and - denotes forwards and backwards spaces respectively.
Definition at line 46 of file library/SolverUtils/RiemannSolvers/UpwindSolver.h.
|
protected |
Default constructor.
Definition at line 66 of file library/SolverUtils/RiemannSolvers/UpwindSolver.cpp.
Referenced by create().
|
inlinestatic |
Definition at line 49 of file library/SolverUtils/RiemannSolvers/UpwindSolver.h.
References UpwindSolver().
|
protectedvirtual |
Implementation of the upwind solver.
The upwind solver assumes that a scalar field Vn is defined, which corresponds with the dot product , where is the advection velocity and defines the normal of a vertex, edge or face at each quadrature point of the trace space.
Fwd | Forwards trace space. |
Bwd | Backwards trace space. |
flux | Resulting flux. |
Implements Nektar::SolverUtils::RiemannSolver.
Definition at line 84 of file library/SolverUtils/RiemannSolvers/UpwindSolver.cpp.
References ASSERTL1, Nektar::SolverUtils::RiemannSolver::CheckScalars(), and Nektar::SolverUtils::RiemannSolver::m_scalars.
|
static |
Definition at line 55 of file library/SolverUtils/RiemannSolvers/UpwindSolver.h.