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

#include <CompressibleSolver.h>

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

Protected Member Functions

 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_ArraySolve (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)
- 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 Attributes

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.
std::map< std::string,
RSScalarFuncType
m_scalars
 Map of scalar function types.
std::map< std::string,
RSVecFuncType
m_vectors
 Map of vector function types.
std::map< std::string,
RSParamFuncType
m_params
 Map of parameter function types.
std::map< std::string,
RSScalarFuncType
m_auxScal
 Map of auxiliary scalar function types.
std::map< std::string,
RSVecFuncType
m_auxVec
 Map of auxiliary vector function types.
Array< OneD, Array< OneD,
NekDouble > > 
m_rotMat
 Rotation matrices for each trace quadrature point.
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_rotStorage
 Rotation storage.

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

Detailed Description

Definition at line 45 of file CompressibleSolver.h.

Constructor & Destructor Documentation

Nektar::CompressibleSolver::CompressibleSolver ( )
protected

Member Function Documentation

virtual void Nektar::CompressibleSolver::v_ArraySolve ( const Array< OneD, const Array< OneD, NekDouble > > &  Fwd,
const Array< OneD, const Array< OneD, NekDouble > > &  Bwd,
Array< OneD, Array< OneD, NekDouble > > &  flux 
)
inlineprotectedvirtual

Reimplemented in Nektar::AverageSolver.

Definition at line 58 of file CompressibleSolver.h.

References ASSERTL0.

Referenced by v_Solve().

{
ASSERTL0(false, "This function should be defined by subclasses.");
}
virtual void Nektar::CompressibleSolver::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 
)
inlineprotectedvirtual

Reimplemented in Nektar::AUSM0Solver, Nektar::AUSM1Solver, Nektar::AUSM2Solver, Nektar::AUSM3Solver, Nektar::ExactSolverToro, Nektar::HLLSolver, Nektar::RoeSolver, Nektar::HLLCSolver, and Nektar::LaxFriedrichsSolver.

Definition at line 66 of file CompressibleSolver.h.

References ASSERTL0.

Referenced by v_Solve().

{
ASSERTL0(false, "This function should be defined by subclasses.");
}
virtual void Nektar::CompressibleSolver::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 
)
inlineprotectedvirtual

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

Definition at line 74 of file CompressibleSolver.h.

References ASSERTL0.

Referenced by v_Solve().

{
ASSERTL0(false, "This function should be defined by subclasses.");
}
void Nektar::CompressibleSolver::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 46 of file CompressibleSolver.cpp.

References m_pointSolve, v_ArraySolve(), v_PointSolve(), and v_PointSolveVisc().

{
{
int expDim = nDim;
int nvariables = Fwd.num_elements();
NekDouble rhouf, rhovf;
// Check if PDE-based SC is used
if (expDim == 1)
{
for (int i = 0; i < Fwd[0].num_elements(); ++i)
{
Fwd [0][i], Fwd [1][i], 0.0, 0.0, Fwd [2][i],
Bwd [0][i], Bwd [1][i], 0.0, 0.0, Bwd [2][i],
flux[0][i], flux[1][i], rhouf, rhovf, flux[2][i]);
}
}
else if (expDim == 2)
{
if (nvariables == expDim+2)
{
for (int i = 0; i < Fwd[0].num_elements(); ++i)
{
Fwd [0][i], Fwd [1][i], Fwd [2][i], 0.0, Fwd [3][i],
Bwd [0][i], Bwd [1][i], Bwd [2][i], 0.0, Bwd [3][i],
flux[0][i], flux[1][i], flux[2][i], rhovf, flux[3][i]);
}
}
if (nvariables > expDim+2)
{
for (int i = 0; i < Fwd[0].num_elements(); ++i)
{
Fwd [0][i], Fwd [1][i], Fwd [2][i], 0.0, Fwd [3][i], Fwd [4][i],
Bwd [0][i], Bwd [1][i], Bwd [2][i], 0.0, Bwd [3][i], Bwd [4][i],
flux[0][i], flux[1][i], flux[2][i], rhovf, flux[3][i], flux[4][i]);
}
}
}
else if (expDim == 3)
{
for (int i = 0; i < Fwd[0].num_elements(); ++i)
{
Fwd [0][i], Fwd [1][i], Fwd [2][i], Fwd [3][i], Fwd [4][i],
Bwd [0][i], Bwd [1][i], Bwd [2][i], Bwd [3][i], Bwd [4][i],
flux[0][i], flux[1][i], flux[2][i], flux[3][i], flux[4][i]);
}
if (nvariables > expDim+2)
{
for (int i = 0; i < Fwd[0].num_elements(); ++i)
{
Fwd [0][i], Fwd [1][i], Fwd [2][i], Fwd [3][i], Fwd [4][i], Fwd [5][i],
Bwd [0][i], Bwd [1][i], Bwd [2][i], Bwd [3][i], Bwd [4][i], Bwd [5][i],
flux[0][i], flux[1][i], flux[2][i], flux[3][i], flux[4][i], flux[5][i]);
}
}
}
}
else
{
v_ArraySolve(Fwd, Bwd, flux);
}
}

Member Data Documentation

bool Nektar::CompressibleSolver::m_pointSolve
protected

Definition at line 48 of file CompressibleSolver.h.

Referenced by Nektar::AverageSolver::AverageSolver(), and v_Solve().