Nektar++
Static Public Member Functions | Static Public Attributes | Protected Member Functions | List of all members
Nektar::AverageSolver Class Reference

#include <AverageSolver.h>

Inheritance diagram for Nektar::AverageSolver:
[legend]

Static Public Member Functions

static RiemannSolverSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession)
 
static RiemannSolverSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession)
 

Static Public Attributes

static std::string solverName
 

Protected Member Functions

 AverageSolver (const LibUtilities::SessionReaderSharedPtr &pSession)
 
virtual void v_ArraySolve (const Array< OneD, const Array< OneD, NekDouble > > &Fwd, const Array< OneD, const Array< OneD, NekDouble > > &Bwd, Array< OneD, Array< OneD, NekDouble > > &flux)
 Average Riemann solver. More...
 
 AverageSolver (const LibUtilities::SessionReaderSharedPtr &pSession)
 
virtual void v_PointSolve (NekDouble hL, NekDouble huL, NekDouble hvL, NekDouble hR, NekDouble huR, NekDouble hvR, NekDouble &hf, NekDouble &huf, NekDouble &hvf)
 Average Value Riemann solver for the Nonlinear Shallow Water Equations. More...
 
- Protected Member Functions inherited from Nektar::NonlinearSWESolver
 NonlinearSWESolver (const LibUtilities::SessionReaderSharedPtr &pSession)
 
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)
 
- Protected Member Functions inherited from Nektar::SolverUtils::RiemannSolver
SOLVER_UTILS_EXPORT RiemannSolver (const LibUtilities::SessionReaderSharedPtr &pSession)
 
virtual 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...
 
- Protected Member Functions inherited from Nektar::CompressibleSolver
 CompressibleSolver (const LibUtilities::SessionReaderSharedPtr &pSession)
 
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 rhoL, NekDouble rhouL, NekDouble rhovL, NekDouble rhowL, NekDouble EL, NekDouble rhoR, NekDouble rhouR, NekDouble rhovR, NekDouble rhowR, NekDouble ER, NekDouble &rhof, NekDouble &rhouf, NekDouble &rhovf, NekDouble &rhowf, NekDouble &Ef)
 
virtual void v_PointSolveVisc (NekDouble rhoL, NekDouble rhouL, NekDouble rhovL, NekDouble rhowL, NekDouble EL, NekDouble EpsL, NekDouble rhoR, NekDouble rhouR, NekDouble rhovR, NekDouble rhowR, NekDouble ER, NekDouble EpsR, NekDouble &rhof, NekDouble &rhouf, NekDouble &rhovf, NekDouble &rhowf, NekDouble &Ef, NekDouble &Epsf)
 
NekDouble GetRoeSoundSpeed (NekDouble rhoL, NekDouble pL, NekDouble eL, NekDouble HL, NekDouble srL, NekDouble rhoR, NekDouble pR, NekDouble eR, NekDouble HR, NekDouble srR, NekDouble HRoe, NekDouble URoe2, NekDouble srLR)
 

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::NonlinearSWESolver
bool m_pointSolve
 
- 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, RSScalarFuncTypem_scalars
 Map of scalar function types. More...
 
std::map< std::string, RSVecFuncTypem_vectors
 Map of vector function types. More...
 
std::map< std::string, RSParamFuncTypem_params
 Map of parameter function types. More...
 
std::map< std::string, RSScalarFuncTypem_auxScal
 Map of auxiliary scalar function types. More...
 
std::map< std::string, RSVecFuncTypem_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...
 
- Protected Attributes inherited from Nektar::CompressibleSolver
bool m_pointSolve
 
EquationOfStateSharedPtr m_eos
 
bool m_idealGas
 

Detailed Description

Definition at line 42 of file CompressibleFlowSolver/RiemannSolvers/AverageSolver.h.

Constructor & Destructor Documentation

◆ AverageSolver() [1/2]

Nektar::AverageSolver::AverageSolver ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Definition at line 47 of file CompressibleFlowSolver/RiemannSolvers/AverageSolver.cpp.

References Nektar::CompressibleSolver::m_pointSolve.

Referenced by create().

49  : CompressibleSolver(pSession)
50  {
51  m_pointSolve = false;
52  }
CompressibleSolver(const LibUtilities::SessionReaderSharedPtr &pSession)

◆ AverageSolver() [2/2]

Nektar::AverageSolver::AverageSolver ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Member Function Documentation

◆ create() [1/2]

static RiemannSolverSharedPtr Nektar::AverageSolver::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

Definition at line 45 of file CompressibleFlowSolver/RiemannSolvers/AverageSolver.h.

References AverageSolver().

47  {
49  new AverageSolver(pSession));
50  }
std::shared_ptr< RiemannSolver > RiemannSolverSharedPtr
A shared pointer to an EquationSystem object.
AverageSolver(const LibUtilities::SessionReaderSharedPtr &pSession)

◆ create() [2/2]

static RiemannSolverSharedPtr Nektar::AverageSolver::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

Definition at line 47 of file ShallowWaterSolver/RiemannSolvers/AverageSolver.h.

References AverageSolver(), solverName, and v_PointSolve().

49  {
51  new AverageSolver(pSession));
52  }
std::shared_ptr< RiemannSolver > RiemannSolverSharedPtr
A shared pointer to an EquationSystem object.
AverageSolver(const LibUtilities::SessionReaderSharedPtr &pSession)

◆ v_ArraySolve()

void Nektar::AverageSolver::v_ArraySolve ( const Array< OneD, const Array< OneD, NekDouble > > &  Fwd,
const Array< OneD, const Array< OneD, NekDouble > > &  Bwd,
Array< OneD, Array< OneD, NekDouble > > &  flux 
)
protectedvirtual

Average Riemann solver.

Parameters
rhoLDensity left state.
rhoRDensity right state.
rhouLx-momentum component left state.
rhouRx-momentum component right state.
rhovLy-momentum component left state.
rhovRy-momentum component right state.
rhowLz-momentum component left state.
rhowRz-momentum component right state.
ELEnergy left state.
EREnergy right state.
rhofComputed Riemann flux for density.
rhoufComputed Riemann flux for x-momentum component
rhovfComputed Riemann flux for y-momentum component
rhowfComputed Riemann flux for z-momentum component
EfComputed Riemann flux for energy.

Reimplemented from Nektar::CompressibleSolver.

Definition at line 73 of file CompressibleFlowSolver/RiemannSolvers/AverageSolver.cpp.

References Nektar::StdRegions::eBwd, Nektar::StdRegions::eFwd, and Nektar::CompressibleSolver::m_eos.

77  {
78  int expDim = Fwd.num_elements()-2;
79  int i, j;
80 
81  for (j = 0; j < Fwd[0].num_elements(); ++j)
82  {
83  NekDouble tmp1 = 0.0, tmp2 = 0.0;
84  Array<OneD, NekDouble> Ufwd(expDim);
85  Array<OneD, NekDouble> Ubwd(expDim);
86 
87  for (i = 0; i < expDim; ++i)
88  {
89  Ufwd[i] = Fwd[i+1][j]/Fwd[0][j];
90  Ubwd[i] = Bwd[i+1][j]/Bwd[0][j];
91  tmp1 += Ufwd[i]*Fwd[i+1][j];
92  tmp2 += Ubwd[i]*Bwd[i+1][j];
93  }
94 
95  // Internal energy
96  NekDouble eFwd = (Fwd[expDim+1][j] - 0.5 * tmp1) / Fwd[0][j];
97  NekDouble eBwd = (Bwd[expDim+1][j] - 0.5 * tmp2) / Bwd[0][j];
98  // Pressure
99  NekDouble Pfwd = m_eos->GetPressure(Fwd[0][j], eFwd);
100  NekDouble Pbwd = m_eos->GetPressure(Bwd[0][j], eBwd);
101 
102  // Compute the average flux
103  flux[0][j] = 0.5 * (Fwd[1][j] + Bwd[1][j]);
104  flux[expDim+1][j] = 0.5 * (Ufwd[0] * (Fwd[expDim+1][j] + Pfwd) +
105  Ubwd[0] * (Bwd[expDim+1][j] + Pbwd));
106 
107  for (i = 0; i < expDim; ++i)
108  {
109  flux[i+1][j] = 0.5 * (Fwd[0][j] * Ufwd[0] * Ufwd[i] +
110  Bwd[0][j] * Ubwd[0] * Ubwd[i]);
111  }
112 
113  // Add in pressure contribution to u field
114  flux[1][j] += 0.5 * (Pfwd + Pbwd);
115  }
116  }
double NekDouble
EquationOfStateSharedPtr m_eos

◆ v_PointSolve()

void Nektar::AverageSolver::v_PointSolve ( NekDouble  hL,
NekDouble  huL,
NekDouble  hvL,
NekDouble  hR,
NekDouble  huR,
NekDouble  hvR,
NekDouble hf,
NekDouble huf,
NekDouble hvf 
)
protectedvirtual

Average Value Riemann solver for the Nonlinear Shallow Water Equations.

Parameters
hLWater depth left state.
hRWater depth right state.
huLx-momentum component left state.
huRx-momentum component right state.
hvLy-momentum component left state.
hvRy-momentum component right state.
hfComputed Riemann flux for density.
hufComputed Riemann flux for x-momentum component
hvfComputed Riemann flux for y-momentum component

Reimplemented from Nektar::NonlinearSWESolver.

Definition at line 68 of file ShallowWaterSolver/RiemannSolvers/AverageSolver.cpp.

References Nektar::SolverUtils::RiemannSolver::m_params.

Referenced by create().

72  {
73  static NekDouble g = m_params["gravity"]();
74 
75  // Left and Right velocities
76  NekDouble uL = huL / hL;
77  NekDouble vL = hvL / hL;
78  NekDouble uR = huR / hR;
79  NekDouble vR = hvR / hR;
80 
81  // Compute fluxes
82  hf = 0.5 * (hL * uL + hR * uR);
83  huf = 0.5 * ((uL * uL * hL + 0.5 * g * hL * hL) +
84  (uR * uR * hR + 0.5 * g * hR * hR));
85  hvf = 0.5 * (hL * uL * vL + hR * uR * vR);
86  }
double NekDouble
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.

Member Data Documentation

◆ solverName

static std::string Nektar::AverageSolver::solverName
static
Initial value:

Definition at line 52 of file CompressibleFlowSolver/RiemannSolvers/AverageSolver.h.

Referenced by create().