Nektar++
Protected Member Functions | List of all members
Nektar::AcousticSolver Class Referenceabstract

#include <AcousticSolver.h>

Inheritance diagram for Nektar::AcousticSolver:
[legend]

Protected Member Functions

 AcousticSolver (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 pL, NekDouble rhoL, NekDouble uL, NekDouble vL, NekDouble wL, NekDouble pR, NekDouble rhoR, NekDouble uR, NekDouble vR, NekDouble wR, NekDouble c0sqL, NekDouble rho0L, NekDouble u0L, NekDouble v0L, NekDouble w0L, NekDouble c0sqR, NekDouble rho0R, NekDouble u0R, NekDouble v0R, NekDouble w0R, NekDouble &pF, NekDouble &rhoF, NekDouble &uF, NekDouble &vF, NekDouble &wF)=0
 
void GetRotBasefield (Array< OneD, Array< OneD, NekDouble >> &bfFwd, Array< OneD, Array< OneD, NekDouble >> &bfBwd)
 
- Protected Member Functions inherited from Nektar::SolverUtils::RiemannSolver
SOLVER_UTILS_EXPORT RiemannSolver (const LibUtilities::SessionReaderSharedPtr &pSession)
 
virtual SOLVER_UTILS_EXPORT ~RiemannSolver ()
 
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)=0
 
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, 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...
 

Detailed Description

Definition at line 47 of file AcousticSolver.h.

Constructor & Destructor Documentation

◆ AcousticSolver()

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

Definition at line 46 of file RiemannSolvers/AcousticSolver.cpp.

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

48  : RiemannSolver(pSession)
49 {
50  m_requiresRotation = true;
51 }
SOLVER_UTILS_EXPORT RiemannSolver(const LibUtilities::SessionReaderSharedPtr &pSession)
bool m_requiresRotation
Indicates whether the Riemann solver requires a rotation to be applied to the velocity fields...

Member Function Documentation

◆ GetRotBasefield()

void Nektar::AcousticSolver::GetRotBasefield ( Array< OneD, Array< OneD, NekDouble >> &  bfFwd,
Array< OneD, Array< OneD, NekDouble >> &  bfBwd 
)
protected

Definition at line 117 of file RiemannSolvers/AcousticSolver.cpp.

References ASSERTL1, Nektar::SolverUtils::RiemannSolver::CheckVectors(), Nektar::SolverUtils::RiemannSolver::m_vectors, and Nektar::SolverUtils::RiemannSolver::rotateToNormal().

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

119 {
120  ASSERTL1(CheckVectors("N"), "N not defined.");
121  ASSERTL1(CheckVectors("basefieldFwdBwd"), "basefieldFwdBwd not defined.");
122  const Array<OneD, const Array<OneD, NekDouble>> normals = m_vectors["N"]();
123  const Array<OneD, const Array<OneD, NekDouble>> basefieldFwdBwd =
124  m_vectors["basefieldFwdBwd"]();
125 
126  int nBF = basefieldFwdBwd.num_elements() / 2;
127  int nDim = normals.num_elements();
128 
129  Array<OneD, Array<OneD, NekDouble>> basefieldFwd(nBF);
130  Array<OneD, Array<OneD, NekDouble>> basefieldBwd(nBF);
131 
132  for (int i = 0; i < nBF; i++)
133  {
134  int j = nBF + i;
135  basefieldFwd[i] = basefieldFwdBwd[i];
136  basefieldBwd[i] = basefieldFwdBwd[j];
137  }
138 
139  Array<OneD, Array<OneD, NekDouble>> baseVecLocs(1);
140  baseVecLocs[0] = Array<OneD, NekDouble>(nDim);
141  for (int i = 0; i < nDim; ++i)
142  {
143  baseVecLocs[0][i] = i + 2;
144  }
145  rotateToNormal(basefieldFwd, normals, baseVecLocs, bfFwd);
146  rotateToNormal(basefieldBwd, normals, baseVecLocs, bfBwd);
147 }
SOLVER_UTILS_EXPORT bool CheckVectors(std::string name)
Determine whether a vector has been defined in m_vectors.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
std::map< std::string, RSVecFuncType > m_vectors
Map of vector function types.
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.

◆ v_PointSolve()

virtual void Nektar::AcousticSolver::v_PointSolve ( NekDouble  pL,
NekDouble  rhoL,
NekDouble  uL,
NekDouble  vL,
NekDouble  wL,
NekDouble  pR,
NekDouble  rhoR,
NekDouble  uR,
NekDouble  vR,
NekDouble  wR,
NekDouble  c0sqL,
NekDouble  rho0L,
NekDouble  u0L,
NekDouble  v0L,
NekDouble  w0L,
NekDouble  c0sqR,
NekDouble  rho0R,
NekDouble  u0R,
NekDouble  v0R,
NekDouble  w0R,
NekDouble pF,
NekDouble rhoF,
NekDouble uF,
NekDouble vF,
NekDouble wF 
)
protectedpure virtual

◆ v_Solve()

void Nektar::AcousticSolver::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

Reimplemented in Nektar::LEESolver.

Definition at line 56 of file RiemannSolvers/AcousticSolver.cpp.

References GetRotBasefield(), and v_PointSolve().

60 {
61  int nTracePts = Fwd[0].num_elements();
62 
63  Array<OneD, Array<OneD, NekDouble>> bfFwd(nDim + 2);
64  Array<OneD, Array<OneD, NekDouble>> bfBwd(nDim + 2);
65  for (int i = 0; i < nDim + 2; i++)
66  {
67  bfFwd[i] = Array<OneD, NekDouble>(nTracePts);
68  bfBwd[i] = Array<OneD, NekDouble>(nTracePts);
69  }
70 
71  GetRotBasefield(bfFwd, bfBwd);
72 
73  int expDim = nDim;
74  NekDouble vF, wF, rhoF;
75 
76  if (expDim == 1)
77  {
78  for (int i = 0; i < nTracePts; ++i)
79  {
81  Fwd[0][i], 0.0, Fwd[1][i], 0.0, 0.0,
82  Bwd[0][i], 0.0, Bwd[1][i], 0.0, 0.0,
83  bfFwd[0][i], bfFwd[1][i], bfFwd[2][i], 0.0, 0.0,
84  bfBwd[0][i], bfBwd[1][i], bfBwd[2][i], 0.0, 0.0,
85  flux[0][i], rhoF, flux[1][i], vF, wF);
86  }
87  }
88  else if (expDim == 2)
89  {
90  for (int i = 0; i < nTracePts; ++i)
91  {
93  Fwd[0][i], 0.0, Fwd[1][i], Fwd[2][i], 0.0,
94  Bwd[0][i], 0.0, Bwd[1][i], Bwd[2][i], 0.0,
95  bfFwd[0][i], bfFwd[1][i], bfFwd[2][i], bfFwd[3][i], 0.0,
96  bfBwd[0][i], bfBwd[1][i], bfBwd[2][i], bfBwd[3][i], 0.0,
97  flux[0][i], rhoF, flux[1][i], flux[2][i], wF);
98  }
99  }
100  else if (expDim == 3)
101  {
102  for (int i = 0; i < nTracePts; ++i)
103  {
104  v_PointSolve(
105  Fwd[0][i], 0.0, Fwd[1][i], Fwd[2][i], Fwd[3][i],
106  Bwd[0][i], 0.0, Bwd[1][i], Bwd[2][i], Bwd[3][i],
107  bfFwd[0][i], bfFwd[1][i], bfFwd[2][i], bfFwd[3][i], bfFwd[4][i],
108  bfBwd[0][i], bfBwd[1][i], bfBwd[2][i], bfBwd[3][i], bfBwd[4][i],
109  flux[0][i], rhoF, flux[1][i], flux[2][i], flux[3][i]);
110  }
111  }
112 }
void GetRotBasefield(Array< OneD, Array< OneD, NekDouble >> &bfFwd, Array< OneD, Array< OneD, NekDouble >> &bfBwd)
double NekDouble
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 c0sqL, NekDouble rho0L, NekDouble u0L, NekDouble v0L, NekDouble w0L, NekDouble c0sqR, NekDouble rho0R, NekDouble u0R, NekDouble v0R, NekDouble w0R, NekDouble &pF, NekDouble &rhoF, NekDouble &uF, NekDouble &vF, NekDouble &wF)=0