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::UpwindSolver Class Reference

#include <UpwindSolver.h>

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

Static Public Member Functions

static RiemannSolverSharedPtr create ()
 

Static Public Attributes

static std::string solverName
 

Protected Member Functions

 UpwindSolver ()
 
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)
 Upwind Riemann solver. More...
 
- Protected Member Functions inherited from Nektar::APESolver
 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)
 
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 solvers/APESolver/RiemannSolvers/UpwindSolver.h.

Constructor & Destructor Documentation

Nektar::UpwindSolver::UpwindSolver ( )
protected

Member Function Documentation

static RiemannSolverSharedPtr Nektar::UpwindSolver::create ( )
inlinestatic

Definition at line 51 of file solvers/APESolver/RiemannSolvers/UpwindSolver.h.

52  {
54  }
boost::shared_ptr< RiemannSolver > RiemannSolverSharedPtr
A shared pointer to an EquationSystem object.
void Nektar::UpwindSolver::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 
)
protectedvirtual

Upwind Riemann solver.

Parameters
pLPerturbation pressure left state
rhoLPerturbation density left state
pRPerturbation pressure right state
rhoRPerturbation density right state
uLx perturbation velocity component left state
uRx perturbation velocity component right state
vLy perturbation velocity component left state
vRy perturbation velocity component right state
wLz perturbation velocity component left state
wRz perturbation velocity component right state
p0Base pressure
rho0Base density
u0Base x velocity component
v0Base y velocity component
w0Base z velocity component
pFComputed Riemann flux for perturbation pressure
rhoFComputed Riemann flux for perturbation density
uFComputed Riemann flux for x perturbation velocity component
vFComputed Riemann flux for y perturbation velocity component
wFComputed Riemann flux for z perturbation velocity component

Reimplemented from Nektar::APESolver.

Definition at line 79 of file solvers/APESolver/RiemannSolvers/UpwindSolver.cpp.

References ASSERTL1, Nektar::SolverUtils::RiemannSolver::CheckParams(), and Nektar::SolverUtils::RiemannSolver::m_params.

84 {
85  // fetch params
86  ASSERTL1(CheckParams("Gamma"), "Gamma not defined.");
87  const NekDouble &gamma = m_params["Gamma"]();
88 
89  // Speed of sound
90  NekDouble c = sqrt(gamma * p0 / rho0);
91 
92  Array<OneD, NekDouble> characteristic(4);
93  Array<OneD, NekDouble> W(2);
94  Array<OneD, NekDouble> lambda(2);
95 
96  // compute the wave speeds
97  lambda[0] = u0 + c;
98  lambda[1] = u0 - c;
99 
100  // calculate the caracteristic variables
101  //left characteristics
102  characteristic[0] = pL/2 + uL*c*rho0/2;
103  characteristic[1] = pL/2 - uL*c*rho0/2;
104  //right characteristics
105  characteristic[2] = pR/2 + uR*c*rho0/2;
106  characteristic[3] = pR/2 - uR*c*rho0/2;
107 
108  //take left or right value of characteristic variable
109  for (int j = 0; j < 2; j++)
110  {
111  if (lambda[j] >= 0)
112  {
113  W[j] = characteristic[j];
114  }
115  if (lambda[j] < 0)
116  {
117  W[j] = characteristic[j+2];
118  }
119  }
120 
121  //calculate conservative variables from characteristics
122  NekDouble p = W[0] + W[1];
123  NekDouble u = (W[0] - W[1])/(c*rho0);
124 
125  // assemble the fluxes
126  pF = rho0*u + u0*p/(c*c);
127  uF = p/rho0 + u0*u + v0*vL + w0*wL;
128  vF = 0.0;
129  wF = 0.0;
130 }
double NekDouble
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
SOLVER_UTILS_EXPORT bool CheckParams(std::string name)
Determine whether a parameter has been defined in m_params.
#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::UpwindSolver::solverName
static