Nektar++
Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | List of all members
Nektar::SolverUtils::RiemannSolver Class Referenceabstract

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:
[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. 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)
 
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

int m_spacedim
 

Protected Member Functions

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)
 

Protected Attributes

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

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

Definition at line 57 of file RiemannSolver.h.

Constructor & Destructor Documentation

◆ RiemannSolver() [1/2]

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

Definition at line 80 of file RiemannSolver.cpp.

81{
82}
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_rotStorage
Rotation storage.
bool m_requiresRotation
Indicates whether the Riemann solver requires a rotation to be applied to the velocity fields.

◆ RiemannSolver() [2/2]

Nektar::SolverUtils::RiemannSolver::RiemannSolver ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Definition at line 84 of file RiemannSolver.cpp.

87{
88 boost::ignore_unused(pSession);
89}

◆ ~RiemannSolver()

virtual SOLVER_UTILS_EXPORT Nektar::SolverUtils::RiemannSolver::~RiemannSolver ( )
inlineprotectedvirtual

Definition at line 160 of file RiemannSolver.h.

160{};

Member Function Documentation

◆ CalcFluxJacobian()

void Nektar::SolverUtils::RiemannSolver::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.

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

Definition at line 515 of file RiemannSolver.cpp.

519{
520 int nPts = Fwd[0].size();
521
523 {
524 ASSERTL1(CheckVectors("N"), "N not defined.");
525 ASSERTL1(CheckAuxVec("vecLocs"), "vecLocs not defined.");
526 const Array<OneD, const Array<OneD, NekDouble>> normals =
527 m_vectors["N"]();
528 const Array<OneD, const Array<OneD, NekDouble>> vecLocs =
529 m_auxVec["vecLocs"]();
530
531 v_CalcFluxJacobian(nDim, Fwd, Bwd, normals, FJac, BJac);
532 }
533 else
534 {
535 Array<OneD, Array<OneD, NekDouble>> normals(nDim);
536 for (int i = 0; i < nDim; i++)
537 {
538 normals[i] = Array<OneD, NekDouble>(nPts, 0.0);
539 }
540 Vmath::Fill(nPts, 1.0, normals[0], 1);
541
542 v_CalcFluxJacobian(nDim, Fwd, Bwd, normals, FJac, BJac);
543 }
544}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
SOLVER_UTILS_EXPORT bool CheckAuxVec(std::string name)
Determine whether a vector has been defined in m_auxVec.
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)
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.
std::map< std::string, RSVecFuncType > m_auxVec
Map of auxiliary vector function types.
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:43

References ASSERTL1, CheckAuxVec(), CheckVectors(), Vmath::Fill(), m_auxVec, m_requiresRotation, m_vectors, and v_CalcFluxJacobian().

◆ CheckAuxScal()

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 355 of file RiemannSolver.cpp.

356{
357 return m_auxScal.find(name) != m_auxScal.end();
358}
std::map< std::string, RSScalarFuncType > m_auxScal
Map of auxiliary scalar function types.

References m_auxScal, and CellMLToNektar.pycml::name.

◆ CheckAuxVec()

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 365 of file RiemannSolver.cpp.

366{
367 return m_auxVec.find(name) != m_auxVec.end();
368}

References m_auxVec, and CellMLToNektar.pycml::name.

Referenced by CalcFluxJacobian(), and Solve().

◆ CheckParams()

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 345 of file RiemannSolver.cpp.

346{
347 return m_params.find(name) != m_params.end();
348}
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.

References m_params, and CellMLToNektar.pycml::name.

Referenced by Nektar::UpwindPulseSolver::RiemannSolverUpwind().

◆ CheckScalars()

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 325 of file RiemannSolver.cpp.

326{
327 return m_scalars.find(name) != m_scalars.end();
328}
std::map< std::string, RSScalarFuncType > m_scalars
Map of scalar function types.

References m_scalars, and CellMLToNektar.pycml::name.

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

◆ CheckVectors()

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 335 of file RiemannSolver.cpp.

336{
337 return m_vectors.find(name) != m_vectors.end();
338}

References m_vectors, and CellMLToNektar.pycml::name.

Referenced by CalcFluxJacobian(), Nektar::AcousticSolver::GetRotBasefield(), Solve(), and Nektar::RoeSolverSIMD::v_Solve().

◆ FromToRotation()

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 418 of file RiemannSolver.cpp.

421{
422 NekDouble v[3];
423 NekDouble e, h, f;
424
425 CROSS(v, from, to);
426 e = DOT(from, to);
427 f = (e < 0) ? -e : e;
428 if (f > 1.0 - EPSILON)
429 {
430 NekDouble u[3], v[3];
431 NekDouble x[3];
432 NekDouble c1, c2, c3;
433 int i, j;
434
435 x[0] = (from[0] > 0.0) ? from[0] : -from[0];
436 x[1] = (from[1] > 0.0) ? from[1] : -from[1];
437 x[2] = (from[2] > 0.0) ? from[2] : -from[2];
438
439 if (x[0] < x[1])
440 {
441 if (x[0] < x[2])
442 {
443 x[0] = 1.0;
444 x[1] = x[2] = 0.0;
445 }
446 else
447 {
448 x[2] = 1.0;
449 x[0] = x[1] = 0.0;
450 }
451 }
452 else
453 {
454 if (x[1] < x[2])
455 {
456 x[1] = 1.0;
457 x[0] = x[2] = 0.0;
458 }
459 else
460 {
461 x[2] = 1.0;
462 x[0] = x[1] = 0.0;
463 }
464 }
465
466 u[0] = x[0] - from[0];
467 u[1] = x[1] - from[1];
468 u[2] = x[2] - from[2];
469 v[0] = x[0] - to[0];
470 v[1] = x[1] - to[1];
471 v[2] = x[2] - to[2];
472
473 c1 = 2.0 / DOT(u, u);
474 c2 = 2.0 / DOT(v, v);
475 c3 = c1 * c2 * DOT(u, v);
476
477 for (i = 0; i < 3; i++)
478 {
479 for (j = 0; j < 3; j++)
480 {
481 mat[3 * i + j] =
482 -c1 * u[i] * u[j] - c2 * v[i] * v[j] + c3 * v[i] * u[j];
483 }
484 mat[i + 3 * i] += 1.0;
485 }
486 }
487 else
488 {
489 NekDouble hvx, hvz, hvxy, hvxz, hvyz;
490 h = 1.0 / (1.0 + e);
491 hvx = h * v[0];
492 hvz = h * v[2];
493 hvxy = hvx * v[1];
494 hvxz = hvx * v[2];
495 hvyz = hvz * v[1];
496 mat[0] = e + hvx * v[0];
497 mat[1] = hvxy - v[2];
498 mat[2] = hvxz + v[1];
499 mat[3] = hvxy + v[2];
500 mat[4] = e + h * v[1] * v[1];
501 mat[5] = hvyz - v[0];
502 mat[6] = hvxz - v[1];
503 mat[7] = hvyz + v[0];
504 mat[8] = e + hvz * v[2];
505 }
506}
#define EPSILON
#define CROSS(dest, v1, v2)
#define DOT(v1, v2)
double NekDouble

References CROSS, DOT, and EPSILON.

Referenced by GenerateRotationMatrices().

◆ GenerateRotationMatrices()

void Nektar::SolverUtils::RiemannSolver::GenerateRotationMatrices ( const Array< OneD, const Array< OneD, NekDouble > > &  normals)
protected

Generate rotation matrices for 3D expansions.

Definition at line 373 of file RiemannSolver.cpp.

375{
376 Array<OneD, NekDouble> xdir(3, 0.0);
377 Array<OneD, NekDouble> tn(3);
378 NekDouble tmp[9];
379 const int nq = normals[0].size();
380 int i, j;
381 xdir[0] = 1.0;
382
383 // Allocate storage for rotation matrices.
384 m_rotMat = Array<OneD, Array<OneD, NekDouble>>(9);
385
386 for (i = 0; i < 9; ++i)
387 {
388 m_rotMat[i] = Array<OneD, NekDouble>(nq);
389 }
390 for (i = 0; i < normals[0].size(); ++i)
391 {
392 // Generate matrix which takes us from (1,0,0) vector to trace
393 // normal.
394 tn[0] = normals[0][i];
395 tn[1] = normals[1][i];
396 tn[2] = normals[2][i];
397 FromToRotation(tn, xdir, tmp);
398
399 for (j = 0; j < 9; ++j)
400 {
401 m_rotMat[j][i] = tmp[j];
402 }
403 }
404}
Array< OneD, Array< OneD, NekDouble > > m_rotMat
Rotation matrices for each trace quadrature point.
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.

References FromToRotation(), and m_rotMat.

Referenced by rotateToNormal(), and Nektar::RoeSolverSIMD::v_Solve().

◆ GetParams()

std::map< std::string, RSParamFuncType > & Nektar::SolverUtils::RiemannSolver::GetParams ( )
inline

Definition at line 125 of file RiemannSolver.h.

126 {
127 return m_params;
128 }

References m_params.

◆ GetScalars()

std::map< std::string, RSScalarFuncType > & Nektar::SolverUtils::RiemannSolver::GetScalars ( )
inline

Definition at line 115 of file RiemannSolver.h.

116 {
117 return m_scalars;
118 }

References m_scalars.

◆ GetVectors()

std::map< std::string, RSVecFuncType > & Nektar::SolverUtils::RiemannSolver::GetVectors ( )
inline

Definition at line 120 of file RiemannSolver.h.

121 {
122 return m_vectors;
123 }

References m_vectors.

◆ rotateFromNormal()

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 251 of file RiemannSolver.cpp.

256{
257 for (int i = 0; i < inarray.size(); ++i)
258 {
259 Vmath::Vcopy(inarray[i].size(), inarray[i], 1, outarray[i], 1);
260 }
261
262 for (int i = 0; i < vecLocs.size(); i++)
263 {
264 ASSERTL1(vecLocs[i].size() == normals.size(),
265 "vecLocs[i] element count mismatch");
266
267 switch (normals.size())
268 {
269 case 1:
270 { // do nothing
271 const int nq = normals[0].size();
272 const int vx = (int)vecLocs[i][0];
273 Vmath::Vmul(nq, inarray[vx], 1, normals[0], 1, outarray[vx], 1);
274 break;
275 }
276 case 2:
277 {
278 const int nq = normals[0].size();
279 const int vx = (int)vecLocs[i][0];
280 const int vy = (int)vecLocs[i][1];
281
282 Vmath::Vmul(nq, inarray[vy], 1, normals[1], 1, outarray[vx], 1);
283 Vmath::Vvtvm(nq, inarray[vx], 1, normals[0], 1, outarray[vx], 1,
284 outarray[vx], 1);
285 Vmath::Vmul(nq, inarray[vx], 1, normals[1], 1, outarray[vy], 1);
286 Vmath::Vvtvp(nq, inarray[vy], 1, normals[0], 1, outarray[vy], 1,
287 outarray[vy], 1);
288 break;
289 }
290
291 case 3:
292 {
293 const int nq = normals[0].size();
294 const int vx = (int)vecLocs[i][0];
295 const int vy = (int)vecLocs[i][1];
296 const int vz = (int)vecLocs[i][2];
297
298 Vmath::Vvtvvtp(nq, inarray[vx], 1, m_rotMat[0], 1, inarray[vy],
299 1, m_rotMat[3], 1, outarray[vx], 1);
300 Vmath::Vvtvp(nq, inarray[vz], 1, m_rotMat[6], 1, outarray[vx],
301 1, outarray[vx], 1);
302 Vmath::Vvtvvtp(nq, inarray[vx], 1, m_rotMat[1], 1, inarray[vy],
303 1, m_rotMat[4], 1, outarray[vy], 1);
304 Vmath::Vvtvp(nq, inarray[vz], 1, m_rotMat[7], 1, outarray[vy],
305 1, outarray[vy], 1);
306 Vmath::Vvtvvtp(nq, inarray[vx], 1, m_rotMat[2], 1, inarray[vy],
307 1, m_rotMat[5], 1, outarray[vz], 1);
308 Vmath::Vvtvp(nq, inarray[vz], 1, m_rotMat[8], 1, outarray[vz],
309 1, outarray[vz], 1);
310 break;
311 }
312
313 default:
314 ASSERTL1(false, "Invalid space dimension.");
315 break;
316 }
317 }
318}
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:207
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:569
void Vvtvm(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvm (vector times vector minus vector): z = w*x - y
Definition: Vmath.cpp:593
void Vvtvvtp(int n, const T *v, int incv, const T *w, int incw, const T *x, int incx, const T *y, int incy, T *z, int incz)
vvtvvtp (vector times vector plus vector times vector):
Definition: Vmath.cpp:687
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191

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

Referenced by Solve().

◆ rotateToNormal()

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 167 of file RiemannSolver.cpp.

172{
173 for (int i = 0; i < inarray.size(); ++i)
174 {
175 Vmath::Vcopy(inarray[i].size(), inarray[i], 1, outarray[i], 1);
176 }
177
178 for (int i = 0; i < vecLocs.size(); i++)
179 {
180 ASSERTL1(vecLocs[i].size() == normals.size(),
181 "vecLocs[i] element count mismatch");
182
183 switch (normals.size())
184 {
185 case 1:
186 { // do nothing
187 const int nq = inarray[0].size();
188 const int vx = (int)vecLocs[i][0];
189 Vmath::Vmul(nq, inarray[vx], 1, normals[0], 1, outarray[vx], 1);
190 break;
191 }
192 case 2:
193 {
194 const int nq = inarray[0].size();
195 const int vx = (int)vecLocs[i][0];
196 const int vy = (int)vecLocs[i][1];
197
198 Vmath::Vmul(nq, inarray[vx], 1, normals[0], 1, outarray[vx], 1);
199 Vmath::Vvtvp(nq, inarray[vy], 1, normals[1], 1, outarray[vx], 1,
200 outarray[vx], 1);
201 Vmath::Vmul(nq, inarray[vx], 1, normals[1], 1, outarray[vy], 1);
202 Vmath::Vvtvm(nq, inarray[vy], 1, normals[0], 1, outarray[vy], 1,
203 outarray[vy], 1);
204 break;
205 }
206
207 case 3:
208 {
209 const int nq = inarray[0].size();
210 const int vx = (int)vecLocs[i][0];
211 const int vy = (int)vecLocs[i][1];
212 const int vz = (int)vecLocs[i][2];
213
214 // Generate matrices if they don't already exist.
215 if (m_rotMat.size() == 0)
216 {
218 }
219
220 // Apply rotation matrices.
221 Vmath::Vvtvvtp(nq, inarray[vx], 1, m_rotMat[0], 1, inarray[vy],
222 1, m_rotMat[1], 1, outarray[vx], 1);
223 Vmath::Vvtvp(nq, inarray[vz], 1, m_rotMat[2], 1, outarray[vx],
224 1, outarray[vx], 1);
225 Vmath::Vvtvvtp(nq, inarray[vx], 1, m_rotMat[3], 1, inarray[vy],
226 1, m_rotMat[4], 1, outarray[vy], 1);
227 Vmath::Vvtvp(nq, inarray[vz], 1, m_rotMat[5], 1, outarray[vy],
228 1, outarray[vy], 1);
229 Vmath::Vvtvvtp(nq, inarray[vx], 1, m_rotMat[6], 1, inarray[vy],
230 1, m_rotMat[7], 1, outarray[vz], 1);
231 Vmath::Vvtvp(nq, inarray[vz], 1, m_rotMat[8], 1, outarray[vz],
232 1, outarray[vz], 1);
233 break;
234 }
235
236 default:
237 ASSERTL1(false, "Invalid space dimension.");
238 break;
239 }
240 }
241}
SOLVER_UTILS_EXPORT void GenerateRotationMatrices(const Array< OneD, const Array< OneD, NekDouble > > &normals)
Generate rotation matrices for 3D expansions.

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

Referenced by Nektar::AcousticSolver::GetRotBasefield(), and Solve().

◆ SetAuxScal()

template<typename FuncPointerT , typename ObjectPointerT >
void Nektar::SolverUtils::RiemannSolver::SetAuxScal ( std::string  name,
FuncPointerT  func,
ObjectPointerT  obj 
)
inline

Definition at line 99 of file RiemannSolver.h.

100 {
101 m_auxScal[name] = std::bind(func, obj);
102 }

References m_auxScal, and CellMLToNektar.pycml::name.

◆ SetAuxVec() [1/2]

template<typename FuncPointerT , typename ObjectPointerT >
void Nektar::SolverUtils::RiemannSolver::SetAuxVec ( std::string  name,
FuncPointerT  func,
ObjectPointerT  obj 
)
inline

Definition at line 105 of file RiemannSolver.h.

106 {
107 m_auxVec[name] = std::bind(func, obj);
108 }

References m_auxVec, and CellMLToNektar.pycml::name.

◆ SetAuxVec() [2/2]

void Nektar::SolverUtils::RiemannSolver::SetAuxVec ( std::string  name,
RSVecFuncType  fp 
)
inline

Definition at line 110 of file RiemannSolver.h.

111 {
112 m_auxVec[name] = fp;
113 }

References m_auxVec, and CellMLToNektar.pycml::name.

◆ SetParam() [1/2]

template<typename FuncPointerT , typename ObjectPointerT >
void Nektar::SolverUtils::RiemannSolver::SetParam ( std::string  name,
FuncPointerT  func,
ObjectPointerT  obj 
)
inline

Definition at line 88 of file RiemannSolver.h.

89 {
90 m_params[name] = std::bind(func, obj);
91 }

References m_params, and CellMLToNektar.pycml::name.

◆ SetParam() [2/2]

void Nektar::SolverUtils::RiemannSolver::SetParam ( std::string  name,
RSParamFuncType  fp 
)
inline

Definition at line 93 of file RiemannSolver.h.

94 {
95 m_params[name] = fp;
96 }

References m_params, and CellMLToNektar.pycml::name.

◆ SetScalar() [1/2]

template<typename FuncPointerT , typename ObjectPointerT >
void Nektar::SolverUtils::RiemannSolver::SetScalar ( std::string  name,
FuncPointerT  func,
ObjectPointerT  obj 
)
inline

Definition at line 66 of file RiemannSolver.h.

67 {
68 m_scalars[name] = std::bind(func, obj);
69 }

References m_scalars, and CellMLToNektar.pycml::name.

◆ SetScalar() [2/2]

void Nektar::SolverUtils::RiemannSolver::SetScalar ( std::string  name,
RSScalarFuncType  fp 
)
inline

Definition at line 71 of file RiemannSolver.h.

72 {
73 m_scalars[name] = fp;
74 }

References m_scalars, and CellMLToNektar.pycml::name.

◆ SetVector() [1/2]

template<typename FuncPointerT , typename ObjectPointerT >
void Nektar::SolverUtils::RiemannSolver::SetVector ( std::string  name,
FuncPointerT  func,
ObjectPointerT  obj 
)
inline

Definition at line 77 of file RiemannSolver.h.

78 {
79 m_vectors[name] = std::bind(func, obj);
80 }

References m_vectors, and CellMLToNektar.pycml::name.

◆ SetVector() [2/2]

void Nektar::SolverUtils::RiemannSolver::SetVector ( std::string  name,
RSVecFuncType  fp 
)
inline

Definition at line 82 of file RiemannSolver.h.

83 {
84 m_vectors[name] = fp;
85 }

References m_vectors, and CellMLToNektar.pycml::name.

◆ Solve()

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 107 of file RiemannSolver.cpp.

111{
113 {
114 ASSERTL1(CheckVectors("N"), "N not defined.");
115 ASSERTL1(CheckAuxVec("vecLocs"), "vecLocs not defined.");
116 const Array<OneD, const Array<OneD, NekDouble>> normals =
117 m_vectors["N"]();
118 const Array<OneD, const Array<OneD, NekDouble>> vecLocs =
119 m_auxVec["vecLocs"]();
120
121 int nFields = Fwd.size();
122 int nPts = Fwd[0].size();
123
124 if (m_rotStorage[0].size() != nFields ||
125 m_rotStorage[0][0].size() != nPts)
126 {
127 for (int i = 0; i < 3; ++i)
128 {
129 m_rotStorage[i] = Array<OneD, Array<OneD, NekDouble>>(nFields);
130 for (int j = 0; j < nFields; ++j)
131 {
132 m_rotStorage[i][j] = Array<OneD, NekDouble>(nPts);
133 }
134 }
135 }
136
137 rotateToNormal(Fwd, normals, vecLocs, m_rotStorage[0]);
138 rotateToNormal(Bwd, normals, vecLocs, m_rotStorage[1]);
140 rotateFromNormal(m_rotStorage[2], normals, vecLocs, flux);
141 }
142 else
143 {
144 v_Solve(nDim, Fwd, Bwd, flux);
145 }
146}
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.
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 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.

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

◆ v_CalcFluxJacobian()

void Nektar::SolverUtils::RiemannSolver::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 
)
protectedvirtual

Definition at line 546 of file RiemannSolver.cpp.

551{
552 boost::ignore_unused(nDim, Fwd, Bwd, normals, FJac, BJac);
553 NEKERROR(ErrorUtil::efatal, "v_CalcFluxJacobian not specified.");
554}
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by CalcFluxJacobian().

◆ v_Solve()

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

◆ m_auxScal

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

Map of auxiliary scalar function types.

Definition at line 148 of file RiemannSolver.h.

Referenced by CheckAuxScal(), and SetAuxScal().

◆ m_auxVec

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

Map of auxiliary vector function types.

Definition at line 150 of file RiemannSolver.h.

Referenced by CalcFluxJacobian(), CheckAuxVec(), SetAuxVec(), and Solve().

◆ m_params

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

◆ m_requiresRotation

bool Nektar::SolverUtils::RiemannSolver::m_requiresRotation
protected

◆ m_rotMat

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

Rotation matrices for each trace quadrature point.

Definition at line 152 of file RiemannSolver.h.

Referenced by GenerateRotationMatrices(), rotateFromNormal(), rotateToNormal(), and Nektar::RoeSolverSIMD::v_Solve().

◆ m_rotStorage

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

Rotation storage.

Definition at line 154 of file RiemannSolver.h.

Referenced by Solve().

◆ m_scalars

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

◆ m_spacedim

int Nektar::SolverUtils::RiemannSolver::m_spacedim

Definition at line 130 of file RiemannSolver.h.

◆ m_vectors

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