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

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 59 of file RiemannSolver.h.

Constructor & Destructor Documentation

◆ RiemannSolver()

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

Definition at line 76 of file RiemannSolver.cpp.

78  : m_requiresRotation(false), m_rotStorage (3)
79  {
80  boost::ignore_unused(pSession);
81  }
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()

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

Member Function Documentation

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

References m_auxScal.

Referenced by ~RiemannSolver().

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

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

References m_auxVec.

Referenced by Solve(), and ~RiemannSolver().

375  {
376  return m_auxVec.find(name) != m_auxVec.end();
377  }
std::map< std::string, RSVecFuncType > m_auxVec
Map of auxiliary vector function types.

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

References m_params.

Referenced by Nektar::UpwindPulseSolver::RiemannSolverUpwind(), and ~RiemannSolver().

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

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

References m_scalars.

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

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

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

References m_vectors.

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

345  {
346  return m_vectors.find(name) != m_vectors.end();
347  }
std::map< std::string, RSVecFuncType > m_vectors
Map of vector function types.

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

References CROSS, DOT, and EPSILON.

Referenced by GenerateRotationMatrices(), and ~RiemannSolver().

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

◆ GenerateRotationMatrices()

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

Generate rotation matrices for 3D expansions.

Definition at line 382 of file RiemannSolver.cpp.

References FromToRotation(), and m_rotMat.

Referenced by rotateToNormal(), and ~RiemannSolver().

384  {
385  Array<OneD, NekDouble> xdir(3,0.0);
386  Array<OneD, NekDouble> tn (3);
387  NekDouble tmp[9];
388  const int nq = normals[0].num_elements();
389  int i, j;
390  xdir[0] = 1.0;
391 
392  // Allocate storage for rotation matrices.
393  m_rotMat = Array<OneD, Array<OneD, NekDouble> >(9);
394 
395  for (i = 0; i < 9; ++i)
396  {
397  m_rotMat[i] = Array<OneD, NekDouble>(nq);
398  }
399  for (i = 0; i < normals[0].num_elements(); ++i)
400  {
401  // Generate matrix which takes us from (1,0,0) vector to trace
402  // normal.
403  tn[0] = normals[0][i];
404  tn[1] = normals[1][i];
405  tn[2] = normals[2][i];
406  FromToRotation(tn, xdir, tmp);
407 
408  for (j = 0; j < 9; ++j)
409  {
410  m_rotMat[j][i] = tmp[j];
411  }
412  }
413  }
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...
double NekDouble
Array< OneD, Array< OneD, NekDouble > > m_rotMat
Rotation matrices for each trace quadrature point.

◆ GetParams()

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

Definition at line 133 of file RiemannSolver.h.

References m_params.

134  {
135  return m_params;
136  }
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.

◆ GetScalars()

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

Definition at line 123 of file RiemannSolver.h.

References m_scalars.

124  {
125  return m_scalars;
126  }
std::map< std::string, RSScalarFuncType > m_scalars
Map of scalar function types.

◆ GetVectors()

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

Definition at line 128 of file RiemannSolver.h.

References m_vectors.

129  {
130  return m_vectors;
131  }
std::map< std::string, RSVecFuncType > m_vectors
Map of vector function types.

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

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

Referenced by Solve(), and ~RiemannSolver().

258  {
259  for (int i = 0; i < inarray.num_elements(); ++i)
260  {
261  Vmath::Vcopy(inarray[i].num_elements(), inarray[i], 1,
262  outarray[i], 1);
263  }
264 
265  for (int i = 0; i < vecLocs.num_elements(); i++)
266  {
267  ASSERTL1(vecLocs[i].num_elements() == normals.num_elements(),
268  "vecLocs[i] element count mismatch");
269 
270  switch (normals.num_elements())
271  {
272  case 1:
273  { // do nothing
274  const int nq = normals[0].num_elements();
275  const int vx = (int)vecLocs[i][0];
276  Vmath::Vmul (nq, inarray [vx], 1, normals [0], 1,
277  outarray[vx], 1);
278  break;
279  }
280  case 2:
281  {
282  const int nq = normals[0].num_elements();
283  const int vx = (int)vecLocs[i][0];
284  const int vy = (int)vecLocs[i][1];
285 
286  Vmath::Vmul (nq, inarray [vy], 1, normals [1], 1,
287  outarray[vx], 1);
288  Vmath::Vvtvm(nq, inarray [vx], 1, normals [0], 1,
289  outarray[vx], 1, outarray[vx], 1);
290  Vmath::Vmul (nq, inarray [vx], 1, normals [1], 1,
291  outarray[vy], 1);
292  Vmath::Vvtvp(nq, inarray [vy], 1, normals [0], 1,
293  outarray[vy], 1, outarray[vy], 1);
294  break;
295  }
296 
297  case 3:
298  {
299  const int nq = normals[0].num_elements();
300  const int vx = (int)vecLocs[i][0];
301  const int vy = (int)vecLocs[i][1];
302  const int vz = (int)vecLocs[i][2];
303 
304  Vmath::Vvtvvtp(nq, inarray [vx], 1, m_rotMat[0], 1,
305  inarray [vy], 1, m_rotMat[3], 1,
306  outarray[vx], 1);
307  Vmath::Vvtvp (nq, inarray [vz], 1, m_rotMat[6], 1,
308  outarray[vx], 1, outarray[vx], 1);
309  Vmath::Vvtvvtp(nq, inarray [vx], 1, m_rotMat[1], 1,
310  inarray [vy], 1, m_rotMat[4], 1,
311  outarray[vy], 1);
312  Vmath::Vvtvp (nq, inarray [vz], 1, m_rotMat[7], 1,
313  outarray[vy], 1, outarray[vy], 1);
314  Vmath::Vvtvvtp(nq, inarray [vx], 1, m_rotMat[2], 1,
315  inarray [vy], 1, m_rotMat[5], 1,
316  outarray[vz], 1);
317  Vmath::Vvtvp (nq, inarray [vz], 1, m_rotMat[8], 1,
318  outarray[vz], 1, outarray[vz], 1);
319  break;
320  }
321 
322  default:
323  ASSERTL1(false, "Invalid space dimension.");
324  break;
325  }
326  }
327  }
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:445
Array< OneD, Array< OneD, NekDouble > > m_rotMat
Rotation matrices for each trace quadrature point.
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:540
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 plus vector): z = w*x - y
Definition: Vmath.cpp:468
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
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:186

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

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

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

167  {
168  for (int i = 0; i < inarray.num_elements(); ++i)
169  {
170  Vmath::Vcopy(inarray[i].num_elements(), inarray[i], 1,
171  outarray[i], 1);
172  }
173 
174  for (int i = 0; i < vecLocs.num_elements(); i++)
175  {
176  ASSERTL1(vecLocs[i].num_elements() == normals.num_elements(),
177  "vecLocs[i] element count mismatch");
178 
179  switch (normals.num_elements())
180  {
181  case 1:
182  { // do nothing
183  const int nq = inarray[0].num_elements();
184  const int vx = (int)vecLocs[i][0];
185  Vmath::Vmul (nq, inarray [vx], 1, normals [0], 1,
186  outarray[vx], 1);
187  break;
188  }
189  case 2:
190  {
191  const int nq = inarray[0].num_elements();
192  const int vx = (int)vecLocs[i][0];
193  const int vy = (int)vecLocs[i][1];
194 
195  Vmath::Vmul (nq, inarray [vx], 1, normals [0], 1,
196  outarray[vx], 1);
197  Vmath::Vvtvp(nq, inarray [vy], 1, normals [1], 1,
198  outarray[vx], 1, outarray[vx], 1);
199  Vmath::Vmul (nq, inarray [vx], 1, normals [1], 1,
200  outarray[vy], 1);
201  Vmath::Vvtvm(nq, inarray [vy], 1, normals [0], 1,
202  outarray[vy], 1, outarray[vy], 1);
203  break;
204  }
205 
206  case 3:
207  {
208  const int nq = inarray[0].num_elements();
209  const int vx = (int)vecLocs[i][0];
210  const int vy = (int)vecLocs[i][1];
211  const int vz = (int)vecLocs[i][2];
212 
213  // Generate matrices if they don't already exist.
214  if (m_rotMat.num_elements() == 0)
215  {
216  GenerateRotationMatrices(normals);
217  }
218 
219  // Apply rotation matrices.
220  Vmath::Vvtvvtp(nq, inarray [vx], 1, m_rotMat[0], 1,
221  inarray [vy], 1, m_rotMat[1], 1,
222  outarray[vx], 1);
223  Vmath::Vvtvp (nq, inarray [vz], 1, m_rotMat[2], 1,
224  outarray[vx], 1, outarray[vx], 1);
225  Vmath::Vvtvvtp(nq, inarray [vx], 1, m_rotMat[3], 1,
226  inarray [vy], 1, m_rotMat[4], 1,
227  outarray[vy], 1);
228  Vmath::Vvtvp (nq, inarray [vz], 1, m_rotMat[5], 1,
229  outarray[vy], 1, outarray[vy], 1);
230  Vmath::Vvtvvtp(nq, inarray [vx], 1, m_rotMat[6], 1,
231  inarray [vy], 1, m_rotMat[7], 1,
232  outarray[vz], 1);
233  Vmath::Vvtvp (nq, inarray [vz], 1, m_rotMat[8], 1,
234  outarray[vz], 1, outarray[vz], 1);
235  break;
236  }
237 
238  default:
239  ASSERTL1(false, "Invalid space dimension.");
240  break;
241  }
242  }
243  }
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:445
void GenerateRotationMatrices(const Array< OneD, const Array< OneD, NekDouble > > &normals)
Generate rotation matrices for 3D expansions.
Array< OneD, Array< OneD, NekDouble > > m_rotMat
Rotation matrices for each trace quadrature point.
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:540
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 plus vector): z = w*x - y
Definition: Vmath.cpp:468
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
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:186

◆ SetAuxScal()

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

Definition at line 108 of file RiemannSolver.h.

References m_auxScal, and CellMLToNektar.pycml::name.

111  {
112  m_auxScal[name] = std::bind(func, obj);
113  }
std::map< std::string, RSScalarFuncType > m_auxScal
Map of auxiliary scalar function types.

◆ SetAuxVec()

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

Definition at line 116 of file RiemannSolver.h.

References m_auxVec, and CellMLToNektar.pycml::name.

119  {
120  m_auxVec[name] = std::bind(func, obj);
121  }
std::map< std::string, RSVecFuncType > m_auxVec
Map of auxiliary vector function types.

◆ SetParam() [1/2]

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

Definition at line 95 of file RiemannSolver.h.

References m_params, and CellMLToNektar.pycml::name.

98  {
99  m_params[name] = std::bind(func, obj);
100  }
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.

◆ SetParam() [2/2]

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

Definition at line 102 of file RiemannSolver.h.

References m_params, and CellMLToNektar.pycml::name.

103  {
104  m_params[name] = fp;
105  }
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.

◆ SetScalar() [1/2]

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

Definition at line 69 of file RiemannSolver.h.

References m_scalars, and CellMLToNektar.pycml::name.

72  {
73  m_scalars[name] = std::bind(func, obj);
74  }
std::map< std::string, RSScalarFuncType > m_scalars
Map of scalar function types.

◆ SetScalar() [2/2]

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

Definition at line 76 of file RiemannSolver.h.

References m_scalars, and CellMLToNektar.pycml::name.

77  {
78  m_scalars[name] = fp;
79  }
std::map< std::string, RSScalarFuncType > m_scalars
Map of scalar function types.

◆ SetVector() [1/2]

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

Definition at line 82 of file RiemannSolver.h.

References m_vectors, and CellMLToNektar.pycml::name.

85  {
86  m_vectors[name] = std::bind(func, obj);
87  }
std::map< std::string, RSVecFuncType > m_vectors
Map of vector function types.

◆ SetVector() [2/2]

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

Definition at line 89 of file RiemannSolver.h.

References m_vectors, and CellMLToNektar.pycml::name.

90  {
91  m_vectors[name] = fp;
92  }
std::map< std::string, RSVecFuncType > m_vectors
Map of vector function types.

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

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

104  {
105  if (m_requiresRotation)
106  {
107  ASSERTL1(CheckVectors("N"), "N not defined.");
108  ASSERTL1(CheckAuxVec("vecLocs"), "vecLocs not defined.");
109  const Array<OneD, const Array<OneD, NekDouble> > normals =
110  m_vectors["N"]();
111  const Array<OneD, const Array<OneD, NekDouble> > vecLocs =
112  m_auxVec["vecLocs"]();
113 
114  int nFields = Fwd .num_elements();
115  int nPts = Fwd[0].num_elements();
116 
117  if (m_rotStorage[0].num_elements() != nFields ||
118  m_rotStorage[0][0].num_elements() != nPts)
119  {
120  for (int i = 0; i < 3; ++i)
121  {
122  m_rotStorage[i] =
123  Array<OneD, Array<OneD, NekDouble> >(nFields);
124  for (int j = 0; j < nFields; ++j)
125  {
126  m_rotStorage[i][j] = Array<OneD, NekDouble>(nPts);
127  }
128  }
129  }
130 
131  rotateToNormal (Fwd, normals, vecLocs, m_rotStorage[0]);
132  rotateToNormal (Bwd, normals, vecLocs, m_rotStorage[1]);
133  v_Solve (nDim, m_rotStorage[0], m_rotStorage[1],
134  m_rotStorage[2]);
135  rotateFromNormal(m_rotStorage[2], normals, vecLocs, flux);
136  }
137  else
138  {
139  v_Solve(nDim, Fwd, Bwd, flux);
140  }
141  }
SOLVER_UTILS_EXPORT bool CheckVectors(std::string name)
Determine whether a vector has been defined in m_vectors.
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.
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
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_rotStorage
Rotation storage.
SOLVER_UTILS_EXPORT bool CheckAuxVec(std::string name)
Determine whether a vector has been defined in m_auxVec.
std::map< std::string, RSVecFuncType > m_auxVec
Map of auxiliary vector function types.
bool m_requiresRotation
Indicates whether the Riemann solver requires a rotation to be applied to the velocity fields...
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
std::map< std::string, RSVecFuncType > m_vectors
Map of vector function types.
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.

◆ 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 151 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 153 of file RiemannSolver.h.

Referenced by 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

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

Definition at line 143 of file RiemannSolver.h.

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

◆ m_rotMat

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

Rotation matrices for each trace quadrature point.

Definition at line 155 of file RiemannSolver.h.

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

◆ m_rotStorage

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

Rotation storage.

Definition at line 157 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 138 of file RiemannSolver.h.

◆ m_vectors

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

Map of vector function types.

Definition at line 147 of file RiemannSolver.h.

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