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

#include <AverageSolver.h>

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

Static Public Member Functions

static RiemannSolverSharedPtr create ()
static RiemannSolverSharedPtr create ()

Static Public Attributes

static std::string solverName

Protected Member Functions

 AverageSolver ()
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.
 AverageSolver ()
virtual void v_PointSolve (double hL, double huL, double hvL, double hR, double huR, double hvR, double &hf, double &huf, double &hvf)
 Average Value Riemann solver for the Nonlinear Shallow Water Equations.
- Protected Member Functions inherited from Nektar::NonlinearSWESolver
 NonlinearSWESolver ()
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 ()
void GenerateRotationMatrices (const Array< OneD, const Array< OneD, NekDouble > > &normals)
 Generate rotation matrices for 3D expansions.
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.
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.
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.
bool CheckScalars (std::string name)
 Determine whether a scalar has been defined in m_scalars.
bool CheckVectors (std::string name)
 Determine whether a vector has been defined in m_vectors.
bool CheckParams (std::string name)
 Determine whether a parameter has been defined in m_params.
bool CheckAuxScal (std::string name)
 Determine whether a scalar has been defined in m_auxScal.
bool CheckAuxVec (std::string name)
 Determine whether a vector has been defined in m_auxVec.
- Protected Member Functions inherited from Nektar::CompressibleSolver
 CompressibleSolver ()
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)

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.
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::CompressibleSolver
bool m_pointSolve

Detailed Description

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

Constructor & Destructor Documentation

Nektar::AverageSolver::AverageSolver ( )
protected
Nektar::AverageSolver::AverageSolver ( )
protected

Member Function Documentation

static RiemannSolverSharedPtr Nektar::AverageSolver::create ( )
inlinestatic
static RiemannSolverSharedPtr Nektar::AverageSolver::create ( )
inlinestatic
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 72 of file CompressibleFlowSolver/RiemannSolvers/AverageSolver.cpp.

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

{
static NekDouble gamma = m_params["gamma"]();
int expDim = Fwd.num_elements()-2;
int i, j;
for (j = 0; j < Fwd[0].num_elements(); ++j)
{
NekDouble tmp1 = 0.0, tmp2 = 0.0;
Array<OneD, NekDouble> Ufwd(expDim);
Array<OneD, NekDouble> Ubwd(expDim);
for (i = 0; i < expDim; ++i)
{
Ufwd[i] = Fwd[i+1][j]/Fwd[0][j];
Ubwd[i] = Bwd[i+1][j]/Bwd[0][j];
tmp1 += Ufwd[i]*Fwd[i+1][j];
tmp2 += Ubwd[i]*Bwd[i+1][j];
}
NekDouble Pfwd = (gamma - 1.0) * (Fwd[expDim+1][j] - 0.5 * tmp1);
NekDouble Pbwd = (gamma - 1.0) * (Bwd[expDim+1][j] - 0.5 * tmp2);
// Compute the average flux
flux[0][j] = 0.5 * (Fwd[1][j] + Bwd[1][j]);
flux[expDim+1][j] = 0.5 * (Ufwd[0] * (Fwd[expDim+1][j] + Pfwd) +
Ubwd[0] * (Bwd[expDim+1][j] + Pbwd));
for (i = 0; i < expDim; ++i)
{
flux[i+1][j] = 0.5 * (Fwd[0][j] * Ufwd[0] * Ufwd[i] +
Bwd[0][j] * Ubwd[0] * Ubwd[i]);
}
// Add in pressure contribution to u field
flux[1][j] += 0.5 * (Pfwd + Pbwd);
}
}
void Nektar::AverageSolver::v_PointSolve ( double  hL,
double  huL,
double  hvL,
double  hR,
double  huR,
double  hvR,
double &  hf,
double &  huf,
double &  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 67 of file ShallowWaterSolver/RiemannSolvers/AverageSolver.cpp.

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

{
static NekDouble g = m_params["gravity"]();
// Left and Right velocities
NekDouble uL = huL / hL;
NekDouble vL = hvL / hL;
NekDouble uR = huR / hR;
NekDouble vR = hvR / hR;
// Compute fluxes
hf = 0.5 * (hL * uL + hR * uR);
huf = 0.5 * ((uL * uL * hL + 0.5 * g * hL * hL) +
(uR * uR * hR + 0.5 * g * hR * hR));
hvf = 0.5 * (hL * uL * vL + hR * uR * vR);
}

Member Data Documentation

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

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