Nektar++
Protected Member Functions | Protected Attributes | List of all members
Nektar::LinearSWESolver Class Reference

#include <LinearSWESolver.h>

Inheritance diagram for Nektar::LinearSWESolver:
[legend]

Protected Member Functions

 LinearSWESolver (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) override
 
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 etaL, NekDouble uL, NekDouble vL, NekDouble dL, NekDouble etaR, NekDouble uR, NekDouble vR, NekDouble dR, NekDouble &etaf, NekDouble &uf, NekDouble &vf)
 
- Protected Member Functions inherited from Nektar::SolverUtils::RiemannSolver
SOLVER_UTILS_EXPORT RiemannSolver ()
 
SOLVER_UTILS_EXPORT RiemannSolver (const LibUtilities::SessionReaderSharedPtr &pSession)
 
virtual SOLVER_UTILS_EXPORT ~RiemannSolver ()
 
SOLVER_UTILS_EXPORT 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...
 
virtual SOLVER_UTILS_EXPORT void v_CalcFluxJacobian (const int nDim, const Array< OneD, const Array< OneD, NekDouble >> &Fwd, const Array< OneD, const Array< OneD, NekDouble >> &Bwd, const Array< OneD, const Array< OneD, NekDouble >> &normals, DNekBlkMatSharedPtr &FJac, DNekBlkMatSharedPtr &BJac)
 

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. 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...
 

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)
 
void SetAuxVec (std::string name, RSVecFuncType fp)
 
std::map< std::string, RSScalarFuncType > & GetScalars ()
 
std::map< std::string, RSVecFuncType > & GetVectors ()
 
std::map< std::string, RSParamFuncType > & GetParams ()
 
SOLVER_UTILS_EXPORT void CalcFluxJacobian (const int nDim, const Array< OneD, const Array< OneD, NekDouble >> &Fwd, const Array< OneD, const Array< OneD, NekDouble >> &Bwd, DNekBlkMatSharedPtr &FJac, DNekBlkMatSharedPtr &BJac)
 Calculate the flux jacobian of Fwd and Bwd. More...
 
- Public Attributes inherited from Nektar::SolverUtils::RiemannSolver
int m_spacedim
 

Detailed Description

Definition at line 44 of file LinearSWESolver.h.

Constructor & Destructor Documentation

◆ LinearSWESolver()

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

Definition at line 41 of file LinearSWESolver.cpp.

43  : RiemannSolver(pSession), m_pointSolve(true)
44 {
45  m_requiresRotation = true;
46 }
bool m_requiresRotation
Indicates whether the Riemann solver requires a rotation to be applied to the velocity fields.
SOLVER_UTILS_EXPORT RiemannSolver()

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

Member Function Documentation

◆ v_ArraySolve()

virtual void Nektar::LinearSWESolver::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

Definition at line 56 of file LinearSWESolver.h.

60  {
61  boost::ignore_unused(Fwd, Bwd, flux);
63  "This function should be defined by subclasses.");
64  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by v_Solve().

◆ v_PointSolve()

virtual void Nektar::LinearSWESolver::v_PointSolve ( NekDouble  etaL,
NekDouble  uL,
NekDouble  vL,
NekDouble  dL,
NekDouble  etaR,
NekDouble  uR,
NekDouble  vR,
NekDouble  dR,
NekDouble etaf,
NekDouble uf,
NekDouble vf 
)
inlineprotectedvirtual

Reimplemented in Nektar::LinearHLLSolver, and Nektar::LinearAverageSolver.

Definition at line 66 of file LinearSWESolver.h.

70  {
71  boost::ignore_unused(etaL, uL, vL, dL, etaR, uR, vR, dR, etaf, uf, vf);
73  "This function should be defined by subclasses.");
74  }

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by v_Solve().

◆ v_Solve()

void Nektar::LinearSWESolver::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 
)
overrideprotectedvirtual

Implements Nektar::SolverUtils::RiemannSolver.

Definition at line 48 of file LinearSWESolver.cpp.

52 {
53  boost::ignore_unused(nDim);
54 
55  // extract the forward and backward trace of the depth
56  const Array<OneD, NekDouble> &dFwd = m_scalars["depthFwd"]();
57  const Array<OneD, NekDouble> &dBwd = m_scalars["depthBwd"]();
58 
59  if (m_pointSolve)
60  {
61 
62  int expDim = Fwd.size() - 1;
63  NekDouble vf;
64 
65  if (expDim == 1)
66  {
67  for (int i = 0; i < Fwd[0].size(); ++i)
68  {
69  v_PointSolve(Fwd[0][i], Fwd[1][i], 0.0, dFwd[i], Bwd[0][i],
70  Bwd[1][i], 0.0, dBwd[i], flux[0][i], flux[1][i],
71  vf);
72  }
73  }
74  else if (expDim == 2)
75  {
76  for (int i = 0; i < Fwd[0].size(); ++i)
77  {
78  v_PointSolve(Fwd[0][i], Fwd[1][i], Fwd[2][i], dFwd[i],
79  Bwd[0][i], Bwd[1][i], Bwd[2][i], dBwd[i],
80  flux[0][i], flux[1][i], flux[2][i]);
81  }
82  }
83  else if (expDim == 3)
84  {
85  ASSERTL0(false, "No 3D Shallow water supported.");
86  }
87  }
88  else
89  {
90  v_ArraySolve(Fwd, Bwd, flux);
91  }
92 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
virtual void v_PointSolve(NekDouble etaL, NekDouble uL, NekDouble vL, NekDouble dL, NekDouble etaR, NekDouble uR, NekDouble vR, NekDouble dR, NekDouble &etaf, NekDouble &uf, NekDouble &vf)
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)
std::map< std::string, RSScalarFuncType > m_scalars
Map of scalar function types.
double NekDouble

References ASSERTL0, m_pointSolve, Nektar::SolverUtils::RiemannSolver::m_scalars, v_ArraySolve(), and v_PointSolve().

Member Data Documentation

◆ m_pointSolve

bool Nektar::LinearSWESolver::m_pointSolve
protected

Definition at line 47 of file LinearSWESolver.h.

Referenced by v_Solve().