| 
    Nektar++
    
   | 
 
Upwind scheme Riemann solver. More...
#include <UpwindSolver.h>
Static Public Member Functions | |
| static SOLVER_UTILS_EXPORT RiemannSolverSharedPtr | create (const LibUtilities::SessionReaderSharedPtr &pSession) | 
Static Public Attributes | |
| static std::string | solverName | 
Protected Member Functions | |
| UpwindSolver (const LibUtilities::SessionReaderSharedPtr &pSession) | |
| Default constructor.  More... | |
| 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) override final | 
| Implementation of the upwind solver.  More... | |
  Protected Member Functions inherited from Nektar::SolverUtils::RiemannSolver | |
| SOLVER_UTILS_EXPORT | RiemannSolver () | 
| SOLVER_UTILS_EXPORT | RiemannSolver (const LibUtilities::SessionReaderSharedPtr &pSession) | 
| virtual SOLVER_UTILS_EXPORT | ~RiemannSolver () | 
| SOLVER_UTILS_EXPORT 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... | |
| virtual SOLVER_UTILS_EXPORT void | v_CalcFluxJacobian (const int nDim, const Array< OneD, const Array< OneD, NekDouble >> &Fwd, const Array< OneD, const Array< OneD, NekDouble >> &Bwd, const Array< OneD, const Array< OneD, NekDouble >> &normals, DNekBlkMatSharedPtr &FJac, DNekBlkMatSharedPtr &BJac) | 
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.  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) | 
| void | SetAuxVec (std::string name, RSVecFuncType fp) | 
| std::map< std::string, RSScalarFuncType > & | GetScalars () | 
| std::map< std::string, RSVecFuncType > & | GetVectors () | 
| std::map< std::string, RSParamFuncType > & | GetParams () | 
| SOLVER_UTILS_EXPORT void | CalcFluxJacobian (const int nDim, const Array< OneD, const Array< OneD, NekDouble >> &Fwd, const Array< OneD, const Array< OneD, NekDouble >> &Bwd, DNekBlkMatSharedPtr &FJac, DNekBlkMatSharedPtr &BJac) | 
| Calculate the flux jacobian of Fwd and Bwd.  More... | |
  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.  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... | |
Upwind scheme Riemann solver.
The upwind solver determines the flux based upon an advection velocity \(\mathbf{V}\) and trace normal \(\mathbf{n}\). In particular, the flux for each component of the velocity field is deterined by:
\[ \mathbf{f}(u^+,u^-) = \begin{cases} \mathbf{V}u^+, & \mathbf{V}\cdot\mathbf{n} \geq 0,\\ \mathbf{V}u^-, & \mathbf{V}\cdot\mathbf{n} < 0.\end{cases} \]
Here the superscript + and - denotes forwards and backwards spaces respectively.
Definition at line 45 of file UpwindSolver.h.
      
  | 
  protected | 
Default constructor.
Definition at line 68 of file UpwindSolver.cpp.
Referenced by create().
      
  | 
  inlinestatic | 
Definition at line 48 of file UpwindSolver.h.
References UpwindSolver().
      
  | 
  finaloverrideprotectedvirtual | 
Implementation of the upwind solver.
The upwind solver assumes that a scalar field Vn is defined, which corresponds with the dot product \(\mathbf{V}\cdot\mathbf{n}\), where \(\mathbf{V}\) is the advection velocity and \(\mathbf{n}\) 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 86 of file UpwindSolver.cpp.
References ASSERTL1, Nektar::SolverUtils::RiemannSolver::CheckScalars(), and Nektar::SolverUtils::RiemannSolver::m_scalars.
      
  | 
  static | 
Definition at line 54 of file UpwindSolver.h.