Nektar++
|
The RiemannSolver class provides an abstract interface under which solvers for various Riemann problems can be implemented. More...
#include <RiemannSolver.h>
Public Member Functions | |
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. More... | |
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 | |
int | m_spacedim |
Protected Member Functions | |
SOLVER_UTILS_EXPORT | RiemannSolver (const LibUtilities::SessionReaderSharedPtr &pSession) |
virtual SOLVER_UTILS_EXPORT | ~RiemannSolver () |
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)=0 |
void | GenerateRotationMatrices (const Array< OneD, const Array< OneD, NekDouble > > &normals) |
Generate rotation matrices for 3D expansions. More... | |
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. More... | |
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. More... | |
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. More... | |
SOLVER_UTILS_EXPORT bool | CheckScalars (std::string name) |
Determine whether a scalar has been defined in m_scalars. More... | |
SOLVER_UTILS_EXPORT bool | CheckVectors (std::string name) |
Determine whether a vector has been defined in m_vectors. More... | |
SOLVER_UTILS_EXPORT bool | CheckParams (std::string name) |
Determine whether a parameter has been defined in m_params. More... | |
SOLVER_UTILS_EXPORT bool | CheckAuxScal (std::string name) |
Determine whether a scalar has been defined in m_auxScal. More... | |
SOLVER_UTILS_EXPORT bool | CheckAuxVec (std::string name) |
Determine whether a vector has been defined in m_auxVec. More... | |
Protected Attributes | |
bool | m_requiresRotation |
Indicates whether the Riemann solver requires a rotation to be applied to the velocity fields. More... | |
std::map< std::string, RSScalarFuncType > | m_scalars |
Map of scalar function types. More... | |
std::map< std::string, RSVecFuncType > | m_vectors |
Map of vector function types. More... | |
std::map< std::string, RSParamFuncType > | m_params |
Map of parameter function types. More... | |
std::map< std::string, RSScalarFuncType > | m_auxScal |
Map of auxiliary scalar function types. More... | |
std::map< std::string, RSVecFuncType > | m_auxVec |
Map of auxiliary vector function types. More... | |
Array< OneD, Array< OneD, NekDouble > > | m_rotMat |
Rotation matrices for each trace quadrature point. More... | |
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > | m_rotStorage |
Rotation storage. More... | |
The RiemannSolver class provides an abstract interface under which solvers for various Riemann problems can be implemented.
Definition at line 59 of file RiemannSolver.h.
|
protected |
Definition at line 76 of file RiemannSolver.cpp.
|
inlineprotectedvirtual |
Definition at line 162 of file RiemannSolver.h.
References CheckAuxScal(), CheckAuxVec(), CheckParams(), CheckScalars(), CheckVectors(), FromToRotation(), GenerateRotationMatrices(), CellMLToNektar.pycml::name, rotateFromNormal(), rotateToNormal(), SOLVER_UTILS_EXPORT, and v_Solve().
|
protected |
Determine whether a scalar has been defined in m_auxScal.
name | Scalar name. |
Definition at line 364 of file RiemannSolver.cpp.
References m_auxScal.
Referenced by ~RiemannSolver().
|
protected |
Determine whether a vector has been defined in m_auxVec.
name | Vector name. |
Definition at line 374 of file RiemannSolver.cpp.
References m_auxVec.
Referenced by Solve(), and ~RiemannSolver().
|
protected |
Determine whether a parameter has been defined in m_params.
name | Parameter name. |
Definition at line 354 of file RiemannSolver.cpp.
References m_params.
Referenced by Nektar::UpwindPulseSolver::RiemannSolverUpwind(), and ~RiemannSolver().
|
protected |
Determine whether a scalar has been defined in m_scalars.
name | Scalar name. |
Definition at line 334 of file RiemannSolver.cpp.
References m_scalars.
Referenced by Nektar::UpwindPulseSolver::v_Solve(), Nektar::SolverUtils::UpwindSolver::v_Solve(), Nektar::SolverUtils::UpwindLDGSolver::v_Solve(), and ~RiemannSolver().
|
protected |
Determine whether a vector has been defined in m_vectors.
name | Vector name. |
Definition at line 344 of file RiemannSolver.cpp.
References m_vectors.
Referenced by Nektar::AcousticSolver::GetRotBasefield(), Solve(), and ~RiemannSolver().
|
protected |
A function for creating a rotation matrix that rotates a vector from into another vector to.
Authors: Tomas Möller, John Hughes "Efficiently Building a Matrix to Rotate One Vector to Another" Journal of Graphics Tools, 4(4):1-4, 1999
from | Normalised 3-vector to rotate from. |
to | Normalised 3-vector to rotate to. |
out | Resulting 3x3 rotation matrix (row-major order). |
Definition at line 427 of file RiemannSolver.cpp.
References CROSS, DOT, and EPSILON.
Referenced by GenerateRotationMatrices(), and ~RiemannSolver().
|
protected |
Generate rotation matrices for 3D expansions.
Definition at line 382 of file RiemannSolver.cpp.
References FromToRotation(), and m_rotMat.
Referenced by rotateToNormal(), and ~RiemannSolver().
|
inline |
Definition at line 133 of file RiemannSolver.h.
References m_params.
|
inline |
Definition at line 123 of file RiemannSolver.h.
References m_scalars.
|
inline |
Definition at line 128 of file RiemannSolver.h.
References m_vectors.
|
protected |
Rotate a vector field from trace normal.
This function performs a rotation of the triad of vector components provided in inarray so that the first component aligns with the Cartesian components; it performs the inverse operation of RiemannSolver::rotateToNormal.
Definition at line 253 of file RiemannSolver.cpp.
References ASSERTL1, m_rotMat, Vmath::Vcopy(), Vmath::Vmul(), Vmath::Vvtvm(), Vmath::Vvtvp(), and Vmath::Vvtvvtp().
Referenced by Solve(), and ~RiemannSolver().
|
protected |
Rotate a vector field to trace normal.
This function performs a rotation of a vector so that the first component aligns with the trace normal direction.
The vectors components are stored in inarray. Their locations must be specified in the "vecLocs" array. vecLocs[0] contains the locations of the first vectors components, vecLocs[1] those of the second and so on.
In 2D, this is accomplished through the transform:
\[ (u_x, u_y) = (n_x u_x + n_y u_y, -n_x v_x + n_y v_y) \]
In 3D, we generate a (non-unique) transformation using RiemannSolver::fromToRotation.
Definition at line 162 of file RiemannSolver.cpp.
References ASSERTL1, GenerateRotationMatrices(), m_rotMat, Vmath::Vcopy(), Vmath::Vmul(), Vmath::Vvtvm(), Vmath::Vvtvp(), and Vmath::Vvtvvtp().
Referenced by Nektar::AcousticSolver::GetRotBasefield(), Solve(), and ~RiemannSolver().
|
inline |
Definition at line 108 of file RiemannSolver.h.
References m_auxScal, and CellMLToNektar.pycml::name.
|
inline |
Definition at line 116 of file RiemannSolver.h.
References m_auxVec, and CellMLToNektar.pycml::name.
|
inline |
Definition at line 95 of file RiemannSolver.h.
References m_params, and CellMLToNektar.pycml::name.
|
inline |
Definition at line 102 of file RiemannSolver.h.
References m_params, and CellMLToNektar.pycml::name.
|
inline |
Definition at line 69 of file RiemannSolver.h.
References m_scalars, and CellMLToNektar.pycml::name.
|
inline |
Definition at line 76 of file RiemannSolver.h.
References m_scalars, and CellMLToNektar.pycml::name.
|
inline |
Definition at line 82 of file RiemannSolver.h.
References m_vectors, and CellMLToNektar.pycml::name.
|
inline |
Definition at line 89 of file RiemannSolver.h.
References m_vectors, and CellMLToNektar.pycml::name.
void Nektar::SolverUtils::RiemannSolver::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.
This routine calls the virtual function v_Solve to perform the Riemann solve. If the flag m_requiresRotation is set, then the velocity field is rotated to the normal direction to perform dimensional splitting, and the resulting fluxes are rotated back to the Cartesian directions before being returned. For the Rotation to work, the normal vectors "N" and the location of the vector components in Fwd "vecLocs"must be set via the SetAuxVec() method.
Fwd | Forwards trace space. |
Bwd | Backwards trace space. |
flux | Resultant flux along trace space. |
Definition at line 99 of file RiemannSolver.cpp.
References ASSERTL1, CheckAuxVec(), CheckVectors(), m_auxVec, m_requiresRotation, m_rotStorage, m_vectors, rotateFromNormal(), rotateToNormal(), and v_Solve().
|
protectedpure virtual |
Implemented in Nektar::SolverUtils::UpwindLDGSolver, Nektar::SolverUtils::UpwindSolver, Nektar::CompressibleSolver, Nektar::NonlinearSWESolver, and Nektar::LinearSWESolver.
Referenced by Solve(), and ~RiemannSolver().
|
protected |
Map of auxiliary scalar function types.
Definition at line 151 of file RiemannSolver.h.
Referenced by CheckAuxScal(), and SetAuxScal().
|
protected |
Map of auxiliary vector function types.
Definition at line 153 of file RiemannSolver.h.
Referenced by CheckAuxVec(), SetAuxVec(), and Solve().
|
protected |
Map of parameter function types.
Definition at line 149 of file RiemannSolver.h.
Referenced by CheckParams(), GetParams(), Nektar::CompressibleSolver::GetRoeSoundSpeed(), Nektar::UpwindPulseSolver::RiemannSolverUpwind(), SetParam(), Nektar::HLLCSolver::v_PointSolve(), Nektar::LinearHLLSolver::v_PointSolve(), Nektar::HLLSolver::v_PointSolve(), Nektar::RoeSolver::v_PointSolve(), Nektar::LaxFriedrichsSolver::v_PointSolve(), Nektar::ExactSolverToro::v_PointSolve(), Nektar::LinearAverageSolver::v_PointSolve(), and Nektar::AverageSolver::v_PointSolve().
|
protected |
Indicates whether the Riemann solver requires a rotation to be applied to the velocity fields.
Definition at line 143 of file RiemannSolver.h.
Referenced by Nektar::AcousticSolver::AcousticSolver(), Nektar::CompressibleSolver::CompressibleSolver(), Nektar::LEESolver::LEESolver(), Nektar::LinearSWESolver::LinearSWESolver(), Nektar::NonlinearSWESolver::NonlinearSWESolver(), and Solve().
Rotation matrices for each trace quadrature point.
Definition at line 155 of file RiemannSolver.h.
Referenced by GenerateRotationMatrices(), rotateFromNormal(), and rotateToNormal().
|
protected |
Map of scalar function types.
Definition at line 145 of file RiemannSolver.h.
Referenced by CheckScalars(), GetScalars(), SetScalar(), Nektar::LinearSWESolver::v_Solve(), Nektar::UpwindPulseSolver::v_Solve(), Nektar::SolverUtils::UpwindSolver::v_Solve(), and Nektar::SolverUtils::UpwindLDGSolver::v_Solve().
int Nektar::SolverUtils::RiemannSolver::m_spacedim |
Definition at line 138 of file RiemannSolver.h.
|
protected |
Map of vector function types.
Definition at line 147 of file RiemannSolver.h.
Referenced by CheckVectors(), Nektar::AcousticSolver::GetRotBasefield(), GetVectors(), SetVector(), and Solve().