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

#include <APEUpwindSolver.h>

Inheritance diagram for Nektar::APEUpwindSolver:
Inheritance graph
[legend]
Collaboration diagram for Nektar::APEUpwindSolver:
Collaboration graph
[legend]

Static Public Member Functions

static RiemannSolverSharedPtr create ()

Static Public Attributes

static std::string solverName

Protected Member Functions

 APEUpwindSolver ()
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 Solve1D (const Array< OneD, const Array< OneD, NekDouble > > &Fwd, const Array< OneD, const Array< OneD, NekDouble > > &Bwd, Array< OneD, Array< OneD, NekDouble > > &flux)
- Protected Member Functions inherited from Nektar::SolverUtils::RiemannSolver
SOLVER_UTILS_EXPORT RiemannSolver ()
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

Array< OneD, Array< OneD,
NekDouble > > 
m_rotBasefield
- 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.
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.

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

Detailed Description

Definition at line 47 of file APEUpwindSolver.h.

Constructor & Destructor Documentation

Nektar::APEUpwindSolver::APEUpwindSolver ( )
protected

Definition at line 53 of file APEUpwindSolver.cpp.

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

:
{
// we need out own rotation logic
}

Member Function Documentation

static RiemannSolverSharedPtr Nektar::APEUpwindSolver::create ( )
inlinestatic

Definition at line 50 of file APEUpwindSolver.h.

void Nektar::APEUpwindSolver::Solve1D ( const Array< OneD, const Array< OneD, NekDouble > > &  Fwd,
const Array< OneD, const Array< OneD, NekDouble > > &  Bwd,
Array< OneD, Array< OneD, NekDouble > > &  flux 
)
protected

Definition at line 123 of file APEUpwindSolver.cpp.

References ASSERTL1, Nektar::SolverUtils::RiemannSolver::CheckParams(), Nektar::SolverUtils::RiemannSolver::CheckVectors(), Nektar::SolverUtils::RiemannSolver::m_params, m_rotBasefield, Nektar::SolverUtils::RiemannSolver::m_vectors, Vmath::Smul(), Vmath::Vmul(), Vmath::Vvtvp(), and Vmath::Zero().

Referenced by v_Solve().

{
// fetch params
ASSERTL1(CheckParams("Gamma"), "Gamma not defined.");
ASSERTL1(CheckParams("Rho"), "Rho not defined.");
ASSERTL1(CheckVectors("N"), "N not defined.");
const NekDouble &gamma= m_params["Gamma"]();
const NekDouble &rho = m_params["Rho"]();
const Array<OneD, const Array<OneD, NekDouble> > normals = m_vectors["N"]();
const Array<OneD, const Array<OneD, NekDouble> > basefield = m_rotBasefield;
int dim = normals.num_elements();
int nvar = dim +1;
ASSERTL1(Fwd.num_elements() == nvar, "Fwd malformed.");
ASSERTL1(Bwd.num_elements() == nvar, "Bwd malformed.");
int nTracePts = Fwd[0].num_elements();
Array<OneD, Array<OneD, NekDouble> > upphysfield(2);
for (int i = 0; i < 2; i++)
{
upphysfield[i] = Array<OneD, NekDouble>(nTracePts);
}
for (int i = 0; i < nTracePts; i++)
{
// assign variables
NekDouble pL = Fwd[0][i];
NekDouble uL = Fwd[1][i];
NekDouble pR = Bwd[0][i];
NekDouble uR = Bwd[1][i];
NekDouble p0 = basefield[0][i];
NekDouble u0 = basefield[1][i];
Array<OneD, NekDouble> characteristic(4);
Array<OneD, NekDouble> W(2);
Array<OneD, NekDouble> lambda(2);
// compute the wave speeds
lambda[0]=u0 + sqrt(p0*gamma*rho)/rho;
lambda[1]=u0 - sqrt(p0*gamma*rho)/rho;
// calculate the caracteristic variables
//left characteristics
characteristic[0] = pL/2 + uL*sqrt(p0*gamma*rho)/2;
characteristic[1] = pL/2 - uL*sqrt(p0*gamma*rho)/2;
//right characteristics
characteristic[2] = pR/2 + uR*sqrt(p0*gamma*rho)/2;
characteristic[3] = pR/2 - uR*sqrt(p0*gamma*rho)/2;
//take left or right value of characteristic variable
for (int j = 0; j < 2; j++)
{
if (lambda[j]>=0)
{
W[j]=characteristic[j];
}
if(lambda[j]<0)
{
W[j]=characteristic[j+2];
}
}
//calculate conservative variables from characteristics
upphysfield[0][i] = W[0]+W[1];
// do not divide by zero
if (p0*gamma*rho != 0)
{
upphysfield[1][i] = (W[0]-W[1])/sqrt(p0*gamma*rho);
}
else
{
upphysfield[1][i] = 0.0;
}
}
// compute the fluxes
// flux[0][i] = gamma*p0 * upphysfield[1] + u0*upphysfield[0]
Vmath::Smul(nTracePts, gamma, basefield[0], 1, flux[0], 1);
Vmath::Vmul(nTracePts, upphysfield[1], 1, flux[0], 1, flux[0], 1);
Vmath::Vvtvp(nTracePts, basefield[1], 1, upphysfield[0], 1, flux[0], 1, flux[0], 1);
// flux[1][i] = upphysfield[0]/rho + u0*upphysfield[1];
Vmath::Smul(nTracePts, 1/rho, upphysfield[0], 1, flux[1], 1);
Vmath::Vvtvp(nTracePts, basefield[1], 1, upphysfield[1], 1, flux[1], 1, flux[1], 1);
for (int j = 2; j < nvar; j++)
{
// flux[1][i] = upphysfield[0]/rho + u0*upphysfield[1] + v0 * vL + w0 + wL;
Vmath::Vvtvp(nTracePts, basefield[j], 1, Fwd[j], 1, flux[1], 1, flux[1], 1);
// flux[2/3][i] = 0.0;
Vmath::Zero(nTracePts, flux[j], 1);
}
}
void Nektar::APEUpwindSolver::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

Implements Nektar::SolverUtils::RiemannSolver.

Definition at line 63 of file APEUpwindSolver.cpp.

References ASSERTL1, Nektar::SolverUtils::RiemannSolver::CheckAuxVec(), Nektar::SolverUtils::RiemannSolver::CheckVectors(), Nektar::SolverUtils::RiemannSolver::m_auxVec, m_rotBasefield, Nektar::SolverUtils::RiemannSolver::m_rotStorage, Nektar::SolverUtils::RiemannSolver::m_vectors, Nektar::SolverUtils::RiemannSolver::rotateFromNormal(), Nektar::SolverUtils::RiemannSolver::rotateToNormal(), and Solve1D().

{
ASSERTL1(CheckVectors("N"), "N not defined.");
ASSERTL1(CheckAuxVec ("vecLocs"), "vecLoc not defined.");
ASSERTL1(CheckVectors("basefield"), "basefield not defined.");
const Array<OneD, const Array<OneD, NekDouble> > normals =
m_vectors["N"]();
const Array<OneD, const Array<OneD, NekDouble> > vecLocs =
m_auxVec["vecLocs"]();
const Array<OneD, const Array<OneD, NekDouble> > basefield =
m_vectors["basefield"]();
int nFields = Fwd .num_elements();
int nPts = Fwd[0].num_elements();
// rotate and store basefield
m_rotBasefield = Array<OneD, Array<OneD, NekDouble> > (nDim+1);
for (int i = 0; i < nDim + 1; i++)
{
m_rotBasefield[i] = Array<OneD, NekDouble>(nPts);
}
Array<OneD, Array<OneD, NekDouble> > baseVecLocs(1);
baseVecLocs[0] = Array<OneD, NekDouble>(nDim);
for (int i = 0; i < nDim; ++i)
{
baseVecLocs[0][i] = 1+i;
}
rotateToNormal(basefield, normals, baseVecLocs, m_rotBasefield);
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);
}

Member Data Documentation

Array<OneD, Array<OneD, NekDouble> > Nektar::APEUpwindSolver::m_rotBasefield
protected

Definition at line 71 of file APEUpwindSolver.h.

Referenced by Solve1D(), and v_Solve().

std::string Nektar::APEUpwindSolver::solverName
static
Initial value:
RegisterCreatorFunction("APEUpwind", APEUpwindSolver::create, "APE Upwind solver")

Definition at line 55 of file APEUpwindSolver.h.