Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Static Public Member Functions | Static Public Attributes | Protected Member Functions | List of all members
Nektar::SolverUtils::UpwindLDGSolver Class Reference

Upwind scheme Riemann solver for LDG. More...

#include <UpwindLDGSolver.h>

Inheritance diagram for Nektar::SolverUtils::UpwindLDGSolver:
Inheritance graph
[legend]
Collaboration diagram for Nektar::SolverUtils::UpwindLDGSolver:
Collaboration graph
[legend]

Static Public Member Functions

static SOLVER_UTILS_EXPORT
RiemannSolverSharedPtr 
create ()
 

Static Public Attributes

static std::string solverName
 

Protected Member Functions

 UpwindLDGSolver ()
 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)
 Implementation of the upwind LDG solver. More...
 
- 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. 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...
 

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)
 
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. 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...
 

Detailed Description

Upwind scheme Riemann solver for LDG.

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 46 of file UpwindLDGSolver.h.

Constructor & Destructor Documentation

Nektar::SolverUtils::UpwindLDGSolver::UpwindLDGSolver ( )
protected

Default constructor.

Definition at line 66 of file UpwindLDGSolver.cpp.

Referenced by create().

66  :
68  {
69  }
SOLVER_UTILS_EXPORT RiemannSolver()

Member Function Documentation

static SOLVER_UTILS_EXPORT RiemannSolverSharedPtr Nektar::SolverUtils::UpwindLDGSolver::create ( )
inlinestatic

Definition at line 49 of file UpwindLDGSolver.h.

References UpwindLDGSolver().

50  {
52  new UpwindLDGSolver());
53  }
boost::shared_ptr< RiemannSolver > RiemannSolverSharedPtr
A shared pointer to an EquationSystem object.
UpwindLDGSolver()
Default constructor.
void Nektar::SolverUtils::UpwindLDGSolver::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 
)
protectedvirtual

Implementation of the upwind LDG 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.

Parameters
FwdForwards trace space.
BwdBackwards trace space.
fluxResulting flux.

Implements Nektar::SolverUtils::RiemannSolver.

Definition at line 84 of file UpwindLDGSolver.cpp.

References ASSERTL1, Nektar::SolverUtils::RiemannSolver::CheckScalars(), and Nektar::SolverUtils::RiemannSolver::m_scalars.

89  {
90  ASSERTL1(CheckScalars("Vn"), "Vn not defined.");
91  const Array<OneD, NekDouble> &traceVel = m_scalars["Vn"]();
92 
93  for (int j = 0; j < traceVel.num_elements(); ++j)
94  {
95  const Array<OneD, const Array<OneD, NekDouble> > &tmp =
96  traceVel[j] >= 0 ? Fwd : Bwd;
97  for (int i = 0; i < Fwd.num_elements(); ++i)
98  {
99  flux[i][j] = tmp[i][j];
100  }
101  }
102  }
SOLVER_UTILS_EXPORT bool CheckScalars(std::string name)
Determine whether a scalar has been defined in m_scalars.
std::map< std::string, RSScalarFuncType > m_scalars
Map of scalar function types.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218

Member Data Documentation

std::string Nektar::SolverUtils::UpwindLDGSolver::solverName
static
Initial value:
RegisterCreatorFunction("UpwindLDG", UpwindLDGSolver::create, "Upwind LDG solver")

Definition at line 55 of file UpwindLDGSolver.h.