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)
 
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_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 ()
 
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
 
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)
 

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 SetALEFlag (bool &ALE)
 
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
 
- 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...
 
bool m_ALESolver = false
 Flag if using the ALE formulation. 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.

48 : RiemannSolver(pSession)
49{
50 m_requiresRotation = true;
51}
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

◆ GetRotBasefield()

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

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

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

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

◆ 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 
)
overrideprotectedvirtual

Implements Nektar::SolverUtils::RiemannSolver.

Reimplemented in Nektar::LEESolver.

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

60{
61 int nTracePts = Fwd[0].size();
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 {
80 v_PointSolve(Fwd[0][i], 0.0, Fwd[1][i], 0.0, 0.0, Bwd[0][i], 0.0,
81 Bwd[1][i], 0.0, 0.0, bfFwd[0][i], bfFwd[1][i],
82 bfFwd[2][i], 0.0, 0.0, bfBwd[0][i], bfBwd[1][i],
83 bfBwd[2][i], 0.0, 0.0, flux[0][i], rhoF, flux[1][i],
84 vF, wF);
85 }
86 }
87 else if (expDim == 2)
88 {
89 for (int i = 0; i < nTracePts; ++i)
90 {
91 v_PointSolve(Fwd[0][i], 0.0, Fwd[1][i], Fwd[2][i], 0.0, Bwd[0][i],
92 0.0, Bwd[1][i], Bwd[2][i], 0.0, bfFwd[0][i],
93 bfFwd[1][i], bfFwd[2][i], bfFwd[3][i], 0.0,
94 bfBwd[0][i], bfBwd[1][i], bfBwd[2][i], bfBwd[3][i],
95 0.0, flux[0][i], rhoF, flux[1][i], flux[2][i], wF);
96 }
97 }
98 else if (expDim == 3)
99 {
100 for (int i = 0; i < nTracePts; ++i)
101 {
102 v_PointSolve(Fwd[0][i], 0.0, Fwd[1][i], Fwd[2][i], Fwd[3][i],
103 Bwd[0][i], 0.0, Bwd[1][i], Bwd[2][i], Bwd[3][i],
104 bfFwd[0][i], bfFwd[1][i], bfFwd[2][i], bfFwd[3][i],
105 bfFwd[4][i], bfBwd[0][i], bfBwd[1][i], bfBwd[2][i],
106 bfBwd[3][i], bfBwd[4][i], flux[0][i], rhoF, flux[1][i],
107 flux[2][i], flux[3][i]);
108 }
109 }
110}
void GetRotBasefield(Array< OneD, Array< OneD, NekDouble > > &bfFwd, Array< OneD, Array< OneD, NekDouble > > &bfBwd)
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
double NekDouble

References GetRotBasefield(), and v_PointSolve().