Nektar++
Static Public Member Functions | Static Public Attributes | Protected Member Functions | List of all members
Nektar::UpwindPulseSolver Class Reference

#include <UpwindPulseSolver.h>

Inheritance diagram for Nektar::UpwindPulseSolver:
[legend]

Static Public Member Functions

static RiemannSolverSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession)
 

Static Public Attributes

static std::string solverName
 

Protected Member Functions

 UpwindPulseSolver (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)
 
void RiemannSolverUpwind (NekDouble AL, NekDouble uL, NekDouble AR, NekDouble uR, NekDouble &Aflux, NekDouble &uflux, NekDouble A_0, NekDouble beta, NekDouble n)
 
- 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 44 of file UpwindPulseSolver.h.

Constructor & Destructor Documentation

◆ UpwindPulseSolver()

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

Definition at line 43 of file UpwindPulseSolver.cpp.

45  : RiemannSolver(pSession)
46 {
47 }
SOLVER_UTILS_EXPORT RiemannSolver(const LibUtilities::SessionReaderSharedPtr &pSession)

Member Function Documentation

◆ create()

static RiemannSolverSharedPtr Nektar::UpwindPulseSolver::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

Definition at line 47 of file UpwindPulseSolver.h.

49  {
50  return RiemannSolverSharedPtr(new UpwindPulseSolver(pSession));
51  }
UpwindPulseSolver(const LibUtilities::SessionReaderSharedPtr &pSession)
std::shared_ptr< RiemannSolver > RiemannSolverSharedPtr
A shared pointer to an EquationSystem object.

◆ RiemannSolverUpwind()

void Nektar::UpwindPulseSolver::RiemannSolverUpwind ( NekDouble  AL,
NekDouble  uL,
NekDouble  AR,
NekDouble  uR,
NekDouble Aflux,
NekDouble uflux,
NekDouble  A_0,
NekDouble  beta,
NekDouble  n 
)
protected

Riemann solver for upwinding at an interface between two elements. Uses the characteristic variables for calculating the upwinded state \((A_u,u_u)\) from the left \((A_L,u_L)\) and right state \((A_R,u_R)\). Returns the upwinded flux ${F}^u$ needed for the weak formulation (1). Details can be found in "Pulse wave propagation in the human vascular system", section 3.3

Definition at line 89 of file UpwindPulseSolver.cpp.

References ASSERTL1, Nektar::SolverUtils::RiemannSolver::CheckParams(), Nektar::SolverUtils::RiemannSolver::m_params, and CellMLToNektar.cellml_metadata::p.

Referenced by v_Solve().

94 {
95  Array<OneD, NekDouble> W(2);
96  Array<OneD, NekDouble> upwindedphysfield(2);
97  NekDouble cL = 0.0;
98  NekDouble cR = 0.0;
99  NekDouble p = 0.0;
100  NekDouble p_t = 0.0;
101 
102  ASSERTL1(CheckParams("rho"), "rho not defined.");
103  NekDouble rho = m_params["rho"]();
104 
105  ASSERTL1(CheckParams("pext"), "pext not defined.");
106  NekDouble pext = m_params["pext"]();
107 
108  // Compute the wave speeds. The use of the normal here allows
109  // for the definition of the characteristics to be inverted
110  // (and hence the left and right state) if n is in the -ve
111  // x-direction. This means we end up with the positive
112  // defintion of the flux which has to therefore be multiplied
113  // by the normal at the end of the methods This is a bit of a
114  // mind twister but is efficient from a coding perspective.
115  cL = sqrt(beta * sqrt(AL) / (2 * rho)) * n;
116  cR = sqrt(beta * sqrt(AR) / (2 * rho)) * n;
117 
118  ASSERTL1(fabs(cL + cR) > fabs(uL + uR), "Conditions are not sub-sonic");
119 
120  // If upwinding from left and right for subsonic domain
121  // then know characteristics immediately
122  W[0] = uL + 4 * cL;
123  W[1] = uR - 4 * cR;
124 
125  // Calculate conservative variables from characteristics
126  NekDouble w0mw1 = 0.25 * (W[0] - W[1]);
127  NekDouble fac = rho / (2 * beta);
128  w0mw1 *= w0mw1; // squared
129  w0mw1 *= w0mw1; // fourth power
130  fac *= fac; // squared
131  upwindedphysfield[0] = w0mw1 * fac;
132  upwindedphysfield[1] = 0.5 * (W[0] + W[1]);
133 
134  // Compute the fluxes multipled by the normal.
135  Aflux = upwindedphysfield[0] * upwindedphysfield[1] * n;
136  p = pext + beta * (sqrt(upwindedphysfield[0]) - sqrt(A_0));
137  p_t = 0.5 * (upwindedphysfield[1] * upwindedphysfield[1]) + p / rho;
138  uflux = p_t * n;
139 }
double NekDouble
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
SOLVER_UTILS_EXPORT bool CheckParams(std::string name)
Determine whether a parameter has been defined in m_params.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250

◆ v_Solve()

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

Calculates the third term of the weak form (1): numerical flux at boundary \( \left[ \mathbf{\psi}^{\delta} \cdot \{ \mathbf{F}^u - \mathbf{F}(\mathbf{U}^{\delta}) \} \right]_{x_e^l}^{x_eû} \)

Definition at line 55 of file UpwindPulseSolver.cpp.

References ASSERTL1, Nektar::SolverUtils::RiemannSolver::CheckScalars(), Nektar::SolverUtils::RiemannSolver::m_scalars, and RiemannSolverUpwind().

59 {
60  int i;
61  int nTracePts = Fwd[0].num_elements();
62 
63  ASSERTL1(CheckScalars("A0"), "A0 not defined.");
64  const Array<OneD, NekDouble> &A0 = m_scalars["A0"]();
65 
66  ASSERTL1(CheckScalars("beta"), "beta not defined.");
67  const Array<OneD, NekDouble> &beta = m_scalars["beta"]();
68 
69  ASSERTL1(CheckScalars("N"), "N not defined.");
70  const Array<OneD, NekDouble> &N = m_scalars["N"]();
71 
72  for (i = 0; i < nTracePts; ++i)
73  {
74  RiemannSolverUpwind(Fwd[0][i], Fwd[1][i], Bwd[0][i], Bwd[1][i],
75  flux[0][i], flux[1][i], A0[i], beta[i], N[i]);
76  }
77 }
SOLVER_UTILS_EXPORT bool CheckScalars(std::string name)
Determine whether a scalar has been defined in m_scalars.
std::map< std::string, RSScalarFuncType > m_scalars
Map of scalar function types.
void RiemannSolverUpwind(NekDouble AL, NekDouble uL, NekDouble AR, NekDouble uR, NekDouble &Aflux, NekDouble &uflux, NekDouble A_0, NekDouble beta, NekDouble n)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250

Member Data Documentation

◆ solverName

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

Definition at line 53 of file UpwindPulseSolver.h.