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

The RiemannSolver class provides an abstract interface under which solvers for various Riemann problems can be implemented. More...

#include <RiemannSolver.h>

Inheritance diagram for Nektar::SolverUtils::RiemannSolver:
Inheritance graph
[legend]
Collaboration diagram for Nektar::SolverUtils::RiemannSolver:
Collaboration graph
[legend]

Public Member Functions

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

int m_spacedim

Protected Member Functions

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

Detailed Description

The RiemannSolver class provides an abstract interface under which solvers for various Riemann problems can be implemented.

Definition at line 62 of file RiemannSolver.h.

Constructor & Destructor Documentation

Nektar::SolverUtils::RiemannSolver::RiemannSolver ( )
protected

Definition at line 77 of file RiemannSolver.cpp.

{
}

Member Function Documentation

bool Nektar::SolverUtils::RiemannSolver::CheckAuxScal ( std::string  name)
protected

Determine whether a scalar has been defined in m_auxScal.

Parameters
nameScalar name.

Definition at line 365 of file RiemannSolver.cpp.

References Nektar::iterator, and m_auxScal.

{
m_auxScal.find(name);
return it != m_auxScal.end();
}
bool Nektar::SolverUtils::RiemannSolver::CheckAuxVec ( std::string  name)
protected

Determine whether a vector has been defined in m_auxVec.

Parameters
nameVector name.

Definition at line 378 of file RiemannSolver.cpp.

References Nektar::iterator, and m_auxVec.

Referenced by Solve(), and Nektar::APEUpwindSolver::v_Solve().

{
m_auxVec.find(name);
return it != m_auxVec.end();
}
bool Nektar::SolverUtils::RiemannSolver::CheckParams ( std::string  name)
protected

Determine whether a parameter has been defined in m_params.

Parameters
nameParameter name.

Definition at line 352 of file RiemannSolver.cpp.

References Nektar::iterator, and m_params.

Referenced by Nektar::APEUpwindSolver::Solve1D().

{
m_params.find(name);
return it != m_params.end();
}
bool Nektar::SolverUtils::RiemannSolver::CheckScalars ( std::string  name)
protected

Determine whether a scalar has been defined in m_scalars.

Parameters
nameScalar name.

Definition at line 326 of file RiemannSolver.cpp.

References Nektar::iterator, and m_scalars.

Referenced by Nektar::SolverUtils::UpwindLDGSolver::v_Solve(), and Nektar::SolverUtils::UpwindSolver::v_Solve().

{
m_scalars.find(name);
return it != m_scalars.end();
}
bool Nektar::SolverUtils::RiemannSolver::CheckVectors ( std::string  name)
protected

Determine whether a vector has been defined in m_vectors.

Parameters
nameVector name.

Definition at line 339 of file RiemannSolver.cpp.

References Nektar::iterator, and m_vectors.

Referenced by Solve(), Nektar::APEUpwindSolver::Solve1D(), and Nektar::APEUpwindSolver::v_Solve().

{
m_vectors.find(name);
return it != m_vectors.end();
}
void Nektar::SolverUtils::RiemannSolver::FromToRotation ( Array< OneD, const NekDouble > &  from,
Array< OneD, const NekDouble > &  to,
NekDouble mat 
)
protected

A function for creating a rotation matrix that rotates a vector from into another vector to.

Authors: Tomas Möller, John Hughes "Efficiently Building a Matrix to Rotate One Vector to Another" Journal of Graphics Tools, 4(4):1-4, 1999

Parameters
fromNormalised 3-vector to rotate from.
toNormalised 3-vector to rotate to.
outResulting 3x3 rotation matrix (row-major order).

Definition at line 434 of file RiemannSolver.cpp.

References CROSS, DOT, and EPSILON.

Referenced by GenerateRotationMatrices().

{
NekDouble v[3];
NekDouble e, h, f;
CROSS(v, from, to);
e = DOT(from, to);
f = (e < 0)? -e:e;
if (f > 1.0 - EPSILON)
{
NekDouble u[3], v[3];
NekDouble x[3];
NekDouble c1, c2, c3;
int i, j;
x[0] = (from[0] > 0.0)? from[0] : -from[0];
x[1] = (from[1] > 0.0)? from[1] : -from[1];
x[2] = (from[2] > 0.0)? from[2] : -from[2];
if (x[0] < x[1])
{
if (x[0] < x[2])
{
x[0] = 1.0; x[1] = x[2] = 0.0;
}
else
{
x[2] = 1.0; x[0] = x[1] = 0.0;
}
}
else
{
if (x[1] < x[2])
{
x[1] = 1.0; x[0] = x[2] = 0.0;
}
else
{
x[2] = 1.0; x[0] = x[1] = 0.0;
}
}
u[0] = x[0] - from[0];
u[1] = x[1] - from[1];
u[2] = x[2] - from[2];
v[0] = x[0] - to [0];
v[1] = x[1] - to [1];
v[2] = x[2] - to [2];
c1 = 2.0 / DOT(u, u);
c2 = 2.0 / DOT(v, v);
c3 = c1 * c2 * DOT(u, v);
for (i = 0; i < 3; i++) {
for (j = 0; j < 3; j++) {
mat[3*i+j] = - c1 * u[i] * u[j]
- c2 * v[i] * v[j]
+ c3 * v[i] * u[j];
}
mat[i+3*i] += 1.0;
}
}
else
{
NekDouble hvx, hvz, hvxy, hvxz, hvyz;
h = 1.0/(1.0 + e);
hvx = h * v[0];
hvz = h * v[2];
hvxy = hvx * v[1];
hvxz = hvx * v[2];
hvyz = hvz * v[1];
mat[0] = e + hvx * v[0];
mat[1] = hvxy - v[2];
mat[2] = hvxz + v[1];
mat[3] = hvxy + v[2];
mat[4] = e + h * v[1] * v[1];
mat[5] = hvyz - v[0];
mat[6] = hvxz - v[1];
mat[7] = hvyz + v[0];
mat[8] = e + hvz * v[2];
}
}
void Nektar::SolverUtils::RiemannSolver::GenerateRotationMatrices ( const Array< OneD, const Array< OneD, NekDouble > > &  normals)
protected

Generate rotation matrices for 3D expansions.

Definition at line 389 of file RiemannSolver.cpp.

References FromToRotation(), and m_rotMat.

Referenced by rotateToNormal().

{
Array<OneD, NekDouble> xdir(3,0.0);
Array<OneD, NekDouble> tn (3);
NekDouble tmp[9];
const int nq = normals[0].num_elements();
int i, j;
xdir[0] = 1.0;
// Allocate storage for rotation matrices.
m_rotMat = Array<OneD, Array<OneD, NekDouble> >(9);
for (i = 0; i < 9; ++i)
{
m_rotMat[i] = Array<OneD, NekDouble>(nq);
}
for (i = 0; i < normals[0].num_elements(); ++i)
{
// Generate matrix which takes us from (1,0,0) vector to trace
// normal.
tn[0] = normals[0][i];
tn[1] = normals[1][i];
tn[2] = normals[2][i];
FromToRotation(tn, xdir, tmp);
for (j = 0; j < 9; ++j)
{
m_rotMat[j][i] = tmp[j];
}
}
}
std::map<std::string, RSParamFuncType>& Nektar::SolverUtils::RiemannSolver::GetParams ( )
inline

Definition at line 136 of file RiemannSolver.h.

References m_params.

{
return m_params;
}
std::map<std::string, RSScalarFuncType>& Nektar::SolverUtils::RiemannSolver::GetScalars ( )
inline

Definition at line 126 of file RiemannSolver.h.

References m_scalars.

{
return m_scalars;
}
std::map<std::string, RSVecFuncType>& Nektar::SolverUtils::RiemannSolver::GetVectors ( )
inline

Definition at line 131 of file RiemannSolver.h.

References m_vectors.

{
return m_vectors;
}
void Nektar::SolverUtils::RiemannSolver::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 
)
protected

Rotate a vector field from trace normal.

This function performs a rotation of the triad of vector components provided in inarray so that the first component aligns with the Cartesian components; it performs the inverse operation of RiemannSolver::rotateToNormal.

Definition at line 249 of file RiemannSolver.cpp.

References ASSERTL1, m_rotMat, Vmath::Vcopy(), Vmath::Vmul(), Vmath::Vvtvm(), Vmath::Vvtvp(), and Vmath::Vvtvvtp().

Referenced by Solve(), and Nektar::APEUpwindSolver::v_Solve().

{
for (int i = 0; i < inarray.num_elements(); ++i)
{
Vmath::Vcopy(inarray[i].num_elements(), inarray[i], 1,
outarray[i], 1);
}
for (int i = 0; i < vecLocs.num_elements(); i++)
{
ASSERTL1(vecLocs[i].num_elements() == normals.num_elements(),
"vecLocs[i] element count mismatch");
switch (normals.num_elements())
{
case 1:
// do nothing
break;
case 2:
{
const int nq = normals[0].num_elements();
const int vx = (int)vecLocs[i][0];
const int vy = (int)vecLocs[i][1];
Vmath::Vmul (nq, inarray [vy], 1, normals [1], 1,
outarray[vx], 1);
Vmath::Vvtvm(nq, inarray [vx], 1, normals [0], 1,
outarray[vx], 1, outarray[vx], 1);
Vmath::Vmul (nq, inarray [vx], 1, normals [1], 1,
outarray[vy], 1);
Vmath::Vvtvp(nq, inarray [vy], 1, normals [0], 1,
outarray[vy], 1, outarray[vy], 1);
break;
}
case 3:
{
const int nq = normals[0].num_elements();
const int vx = (int)vecLocs[i][0];
const int vy = (int)vecLocs[i][1];
const int vz = (int)vecLocs[i][2];
Vmath::Vvtvvtp(nq, inarray [vx], 1, m_rotMat[0], 1,
inarray [vy], 1, m_rotMat[3], 1,
outarray[vx], 1);
Vmath::Vvtvp (nq, inarray [vz], 1, m_rotMat[6], 1,
outarray[vx], 1, outarray[vx], 1);
Vmath::Vvtvvtp(nq, inarray [vx], 1, m_rotMat[1], 1,
inarray [vy], 1, m_rotMat[4], 1,
outarray[vy], 1);
Vmath::Vvtvp (nq, inarray [vz], 1, m_rotMat[7], 1,
outarray[vy], 1, outarray[vy], 1);
Vmath::Vvtvvtp(nq, inarray [vx], 1, m_rotMat[2], 1,
inarray [vy], 1, m_rotMat[5], 1,
outarray[vz], 1);
Vmath::Vvtvp (nq, inarray [vz], 1, m_rotMat[8], 1,
outarray[vz], 1, outarray[vz], 1);
break;
}
default:
ASSERTL1(false, "Invalid space dimension.");
break;
}
}
}
void Nektar::SolverUtils::RiemannSolver::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 
)
protected

Rotate a vector field to trace normal.

This function performs a rotation of a vector so that the first component aligns with the trace normal direction.

The vectors components are stored in inarray. Their locations must be specified in the "vecLocs" array. vecLocs[0] contains the locations of the first vectors components, vecLocs[1] those of the second and so on.

In 2D, this is accomplished through the transform:

\[ (u_x, u_y) = (n_x u_x + n_y u_y, -n_x v_x + n_y v_y) \]

In 3D, we generate a (non-unique) transformation using RiemannSolver::fromToRotation.

Definition at line 162 of file RiemannSolver.cpp.

References ASSERTL1, GenerateRotationMatrices(), m_rotMat, Vmath::Vcopy(), Vmath::Vmul(), Vmath::Vvtvm(), Vmath::Vvtvp(), and Vmath::Vvtvvtp().

Referenced by Solve(), and Nektar::APEUpwindSolver::v_Solve().

{
for (int i = 0; i < inarray.num_elements(); ++i)
{
Vmath::Vcopy(inarray[i].num_elements(), inarray[i], 1,
outarray[i], 1);
}
for (int i = 0; i < vecLocs.num_elements(); i++)
{
ASSERTL1(vecLocs[i].num_elements() == normals.num_elements(),
"vecLocs[i] element count mismatch");
switch (normals.num_elements())
{
case 1:
// do nothing
break;
case 2:
{
const int nq = inarray[0].num_elements();
const int vx = (int)vecLocs[i][0];
const int vy = (int)vecLocs[i][1];
Vmath::Vmul (nq, inarray [vx], 1, normals [0], 1,
outarray[vx], 1);
Vmath::Vvtvp(nq, inarray [vy], 1, normals [1], 1,
outarray[vx], 1, outarray[vx], 1);
Vmath::Vmul (nq, inarray [vx], 1, normals [1], 1,
outarray[vy], 1);
Vmath::Vvtvm(nq, inarray [vy], 1, normals [0], 1,
outarray[vy], 1, outarray[vy], 1);
break;
}
case 3:
{
const int nq = inarray[0].num_elements();
const int vx = (int)vecLocs[i][0];
const int vy = (int)vecLocs[i][1];
const int vz = (int)vecLocs[i][2];
// Generate matrices if they don't already exist.
if (m_rotMat.num_elements() == 0)
{
}
// Apply rotation matrices.
Vmath::Vvtvvtp(nq, inarray [vx], 1, m_rotMat[0], 1,
inarray [vy], 1, m_rotMat[1], 1,
outarray[vx], 1);
Vmath::Vvtvp (nq, inarray [vz], 1, m_rotMat[2], 1,
outarray[vx], 1, outarray[vx], 1);
Vmath::Vvtvvtp(nq, inarray [vx], 1, m_rotMat[3], 1,
inarray [vy], 1, m_rotMat[4], 1,
outarray[vy], 1);
Vmath::Vvtvp (nq, inarray [vz], 1, m_rotMat[5], 1,
outarray[vy], 1, outarray[vy], 1);
Vmath::Vvtvvtp(nq, inarray [vx], 1, m_rotMat[6], 1,
inarray [vy], 1, m_rotMat[7], 1,
outarray[vz], 1);
Vmath::Vvtvp (nq, inarray [vz], 1, m_rotMat[8], 1,
outarray[vz], 1, outarray[vz], 1);
break;
}
default:
ASSERTL1(false, "Invalid space dimension.");
break;
}
}
}
template<typename FuncPointerT , typename ObjectPointerT >
void Nektar::SolverUtils::RiemannSolver::SetAuxScal ( std::string  name,
FuncPointerT  func,
ObjectPointerT  obj 
)
inline

Definition at line 111 of file RiemannSolver.h.

References m_auxScal.

{
m_auxScal[name] = boost::bind(func, obj);
}
template<typename FuncPointerT , typename ObjectPointerT >
void Nektar::SolverUtils::RiemannSolver::SetAuxVec ( std::string  name,
FuncPointerT  func,
ObjectPointerT  obj 
)
inline

Definition at line 119 of file RiemannSolver.h.

References m_auxVec.

{
m_auxVec[name] = boost::bind(func, obj);
}
template<typename FuncPointerT , typename ObjectPointerT >
void Nektar::SolverUtils::RiemannSolver::SetParam ( std::string  name,
FuncPointerT  func,
ObjectPointerT  obj 
)
inline

Definition at line 98 of file RiemannSolver.h.

References m_params.

{
m_params[name] = boost::bind(func, obj);
}
void Nektar::SolverUtils::RiemannSolver::SetParam ( std::string  name,
RSParamFuncType  fp 
)
inline

Definition at line 105 of file RiemannSolver.h.

References m_params.

{
m_params[name] = fp;
}
template<typename FuncPointerT , typename ObjectPointerT >
void Nektar::SolverUtils::RiemannSolver::SetScalar ( std::string  name,
FuncPointerT  func,
ObjectPointerT  obj 
)
inline

Definition at line 72 of file RiemannSolver.h.

References m_scalars.

{
m_scalars[name] = boost::bind(func, obj);
}
void Nektar::SolverUtils::RiemannSolver::SetScalar ( std::string  name,
RSScalarFuncType  fp 
)
inline

Definition at line 79 of file RiemannSolver.h.

References m_scalars.

{
m_scalars[name] = fp;
}
template<typename FuncPointerT , typename ObjectPointerT >
void Nektar::SolverUtils::RiemannSolver::SetVector ( std::string  name,
FuncPointerT  func,
ObjectPointerT  obj 
)
inline

Definition at line 85 of file RiemannSolver.h.

References m_vectors.

{
m_vectors[name] = boost::bind(func, obj);
}
void Nektar::SolverUtils::RiemannSolver::SetVector ( std::string  name,
RSVecFuncType  fp 
)
inline

Definition at line 92 of file RiemannSolver.h.

References m_vectors.

{
m_vectors[name] = fp;
}
void Nektar::SolverUtils::RiemannSolver::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.

This routine calls the virtual function v_Solve to perform the Riemann solve. If the flag m_requiresRotation is set, then the velocity field is rotated to the normal direction to perform dimensional splitting, and the resulting fluxes are rotated back to the Cartesian directions before being returned. For the Rotation to work, the normal vectors "N" and the location of the vector components in Fwd "vecLocs"must be set via the SetAuxVec() method.

Parameters
FwdForwards trace space.
BwdBackwards trace space.
fluxResultant flux along trace space.

Definition at line 99 of file RiemannSolver.cpp.

References ASSERTL1, CheckAuxVec(), CheckVectors(), m_auxVec, m_requiresRotation, m_rotStorage, m_vectors, rotateFromNormal(), rotateToNormal(), and v_Solve().

{
{
ASSERTL1(CheckVectors("N"), "N not defined.");
ASSERTL1(CheckAuxVec("vecLocs"), "vecLocs not defined.");
const Array<OneD, const Array<OneD, NekDouble> > normals =
m_vectors["N"]();
const Array<OneD, const Array<OneD, NekDouble> > vecLocs =
m_auxVec["vecLocs"]();
int nFields = Fwd .num_elements();
int nPts = Fwd[0].num_elements();
if (m_rotStorage[0].num_elements() != nFields ||
m_rotStorage[0][0].num_elements() != nPts)
{
for (int i = 0; i < 3; ++i)
{
Array<OneD, Array<OneD, NekDouble> >(nFields);
for (int j = 0; j < nFields; ++j)
{
m_rotStorage[i][j] = Array<OneD, NekDouble>(nPts);
}
}
}
rotateToNormal (Fwd, normals, vecLocs, m_rotStorage[0]);
rotateToNormal (Bwd, normals, vecLocs, m_rotStorage[1]);
rotateFromNormal(m_rotStorage[2], normals, vecLocs, flux);
}
else
{
v_Solve(nDim, Fwd, Bwd, flux);
}
}
virtual void Nektar::SolverUtils::RiemannSolver::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 
)
protectedpure virtual

Member Data Documentation

std::map<std::string, RSScalarFuncType> Nektar::SolverUtils::RiemannSolver::m_auxScal
protected

Map of auxiliary scalar function types.

Definition at line 154 of file RiemannSolver.h.

Referenced by CheckAuxScal(), and SetAuxScal().

std::map<std::string, RSVecFuncType> Nektar::SolverUtils::RiemannSolver::m_auxVec
protected

Map of auxiliary vector function types.

Definition at line 156 of file RiemannSolver.h.

Referenced by CheckAuxVec(), SetAuxVec(), Solve(), and Nektar::APEUpwindSolver::v_Solve().

std::map<std::string, RSParamFuncType > Nektar::SolverUtils::RiemannSolver::m_params
protected
bool Nektar::SolverUtils::RiemannSolver::m_requiresRotation
protected

Indicates whether the Riemann solver requires a rotation to be applied to the velocity fields.

Definition at line 146 of file RiemannSolver.h.

Referenced by Nektar::APEUpwindSolver::APEUpwindSolver(), Nektar::CompressibleSolver::CompressibleSolver(), Nektar::LinearSWESolver::LinearSWESolver(), Nektar::NonlinearSWESolver::NonlinearSWESolver(), and Solve().

Array<OneD, Array<OneD, NekDouble> > Nektar::SolverUtils::RiemannSolver::m_rotMat
protected

Rotation matrices for each trace quadrature point.

Definition at line 158 of file RiemannSolver.h.

Referenced by GenerateRotationMatrices(), rotateFromNormal(), and rotateToNormal().

Array<OneD, Array<OneD, Array<OneD, NekDouble> > > Nektar::SolverUtils::RiemannSolver::m_rotStorage
protected

Rotation storage.

Definition at line 160 of file RiemannSolver.h.

Referenced by Solve(), and Nektar::APEUpwindSolver::v_Solve().

std::map<std::string, RSScalarFuncType> Nektar::SolverUtils::RiemannSolver::m_scalars
protected
int Nektar::SolverUtils::RiemannSolver::m_spacedim

Definition at line 141 of file RiemannSolver.h.

std::map<std::string, RSVecFuncType> Nektar::SolverUtils::RiemannSolver::m_vectors
protected

Map of vector function types.

Definition at line 150 of file RiemannSolver.h.

Referenced by CheckVectors(), GetVectors(), SetVector(), Solve(), Nektar::APEUpwindSolver::Solve1D(), and Nektar::APEUpwindSolver::v_Solve().