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

#include <APESolver.h>

Inheritance diagram for Nektar::APESolver:
Inheritance graph
[legend]
Collaboration diagram for Nektar::APESolver:
Collaboration graph
[legend]

Protected Member Functions

 APESolver ()
 
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)
 
virtual void v_PointSolve (NekDouble pL, NekDouble rhoL, NekDouble uL, NekDouble vL, NekDouble wL, NekDouble pR, NekDouble rhoR, NekDouble uR, NekDouble vR, NekDouble wR, NekDouble p0, NekDouble rho0, NekDouble u0, NekDouble v0, NekDouble w0, NekDouble &pF, NekDouble &rhoF, NekDouble &uF, NekDouble &vF, NekDouble &wF)
 
Array< OneD, Array< OneD,
NekDouble > > 
GetRotBasefield ()
 
- 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

Definition at line 48 of file APESolver.h.

Constructor & Destructor Documentation

Nektar::APESolver::APESolver ( )
protected

Definition at line 47 of file RiemannSolvers/APESolver.cpp.

References Nektar::SolverUtils::RiemannSolver::m_requiresRotation.

47  :
49 {
50  m_requiresRotation = true;
51 }
SOLVER_UTILS_EXPORT RiemannSolver()
bool m_requiresRotation
Indicates whether the Riemann solver requires a rotation to be applied to the velocity fields...

Member Function Documentation

Array< OneD, Array< OneD, NekDouble > > Nektar::APESolver::GetRotBasefield ( )
protected

Definition at line 107 of file RiemannSolvers/APESolver.cpp.

References ASSERTL1, Nektar::SolverUtils::RiemannSolver::CheckVectors(), Nektar::SolverUtils::RiemannSolver::m_vectors, and Nektar::SolverUtils::RiemannSolver::rotateToNormal().

Referenced by v_Solve().

108 {
109  ASSERTL1(CheckVectors("N"), "N not defined.");
110  ASSERTL1(CheckVectors("basefield"), "basefield not defined.");
111  const Array<OneD, const Array<OneD, NekDouble> > normals = m_vectors["N"]();
112  const Array<OneD, const Array<OneD, NekDouble> > basefield =
113  m_vectors["basefield"]();
114 
115  int nTracePts = normals[0].num_elements();
116  int nDim = normals.num_elements();
117 
118  Array< OneD, Array< OneD, NekDouble > > rotBasefield(nDim+2);
119  for (int i = 0; i < nDim + 2; i++)
120  {
121  rotBasefield[i] = Array<OneD, NekDouble>(nTracePts);
122  }
123  Array<OneD, Array<OneD, NekDouble> > baseVecLocs(1);
124  baseVecLocs[0] = Array<OneD, NekDouble>(nDim);
125  for (int i = 0; i < nDim; ++i)
126  {
127  baseVecLocs[0][i] = i + 2;
128  }
129  rotateToNormal(basefield, normals, baseVecLocs, rotBasefield);
130 
131  return rotBasefield;
132 }
SOLVER_UTILS_EXPORT bool CheckVectors(std::string name)
Determine whether a vector has been defined in m_vectors.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
std::map< std::string, RSVecFuncType > m_vectors
Map of vector function types.
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.
virtual void Nektar::APESolver::v_PointSolve ( NekDouble  pL,
NekDouble  rhoL,
NekDouble  uL,
NekDouble  vL,
NekDouble  wL,
NekDouble  pR,
NekDouble  rhoR,
NekDouble  uR,
NekDouble  vR,
NekDouble  wR,
NekDouble  p0,
NekDouble  rho0,
NekDouble  u0,
NekDouble  v0,
NekDouble  w0,
NekDouble pF,
NekDouble rhoF,
NekDouble uF,
NekDouble vF,
NekDouble wF 
)
inlineprotectedvirtual

Reimplemented in Nektar::LaxFriedrichsSolver, and Nektar::UpwindSolver.

Definition at line 59 of file APESolver.h.

References ASSERTL0.

Referenced by v_Solve().

64  {
65  ASSERTL0(false, "This function should be defined by subclasses.");
66  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
void Nektar::APESolver::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

Implements Nektar::SolverUtils::RiemannSolver.

Definition at line 57 of file RiemannSolvers/APESolver.cpp.

References GetRotBasefield(), and v_PointSolve().

62 {
63  Array< OneD, Array< OneD, NekDouble > > bf = GetRotBasefield();
64 
65  int expDim = nDim;
66  NekDouble vF, wF, rhoF;
67 
68  if (expDim == 1)
69  {
70  for (int i = 0; i < Fwd[0].num_elements(); ++i)
71  {
73  Fwd[0][i], 0.0, Fwd[1][i], 0.0, 0.0,
74  Bwd[0][i], 0.0, Bwd[1][i], 0.0, 0.0,
75  bf[0][i], bf[1][i], bf[2][i], 0.0, 0.0,
76  flux[0][i], rhoF, flux[1][i], vF, wF);
77  }
78  }
79  else if (expDim == 2)
80  {
81  for (int i = 0; i < Fwd[0].num_elements(); ++i)
82  {
84  Fwd[0][i], 0.0, Fwd[1][i], Fwd[2][i], 0.0,
85  Bwd[0][i], 0.0, Bwd[1][i], Bwd[2][i], 0.0,
86  bf[0][i], bf[1][i], bf[2][i], bf[3][i], 0.0,
87  flux[0][i], rhoF, flux[1][i], flux[2][i], wF);
88  }
89  }
90  else if (expDim == 3)
91  {
92  for (int i = 0; i < Fwd[0].num_elements(); ++i)
93  {
95  Fwd[0][i], 0.0, Fwd[1][i], Fwd[2][i], Fwd[3][i],
96  Bwd[0][i], 0.0, Bwd[1][i], Bwd[2][i], Bwd[3][i],
97  bf[0][i], bf[1][i], bf[2][i], bf[3][i], bf[4][i],
98  flux[0][i], rhoF, flux[1][i], flux[2][i], flux[3][i]);
99  }
100  }
101 }
Array< OneD, Array< OneD, NekDouble > > GetRotBasefield()
double NekDouble
virtual void v_PointSolve(NekDouble pL, NekDouble rhoL, NekDouble uL, NekDouble vL, NekDouble wL, NekDouble pR, NekDouble rhoR, NekDouble uR, NekDouble vR, NekDouble wR, NekDouble p0, NekDouble rho0, NekDouble u0, NekDouble v0, NekDouble w0, NekDouble &pF, NekDouble &rhoF, NekDouble &uF, NekDouble &vF, NekDouble &wF)
Definition: APESolver.h:59