Nektar++
RiemannSolver.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: RiemannSolver.cpp
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: Abstract base class for Riemann solvers with factory.
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <boost/core/ignore_unused.hpp>
36 
39 
40 #define EPSILON 0.000001
41 
42 #define CROSS(dest, v1, v2){ \
43  dest[0] = v1[1] * v2[2] - v1[2] * v2[1]; \
44  dest[1] = v1[2] * v2[0] - v1[0] * v2[2]; \
45  dest[2] = v1[0] * v2[1] - v1[1] * v2[0];}
46 
47 #define DOT(v1, v2) (v1[0] * v2[0] + v1[1] * v2[1] + v1[2] * v2[2])
48 
49 #define SUB(dest, v1, v2){ \
50  dest[0] = v1[0] - v2[0]; \
51  dest[1] = v1[1] - v2[1]; \
52  dest[2] = v1[2] - v2[2];}
53 
54 using namespace std;
55 
56 namespace Nektar
57 {
58  namespace SolverUtils
59  {
60  /**
61  * Retrieves the singleton instance of the Riemann solver factory.
62  */
64  {
65  static RiemannSolverFactory instance;
66  return instance;
67  }
68 
69  /**
70  * @class RiemannSolver
71  *
72  * @brief The RiemannSolver class provides an abstract interface under
73  * which solvers for various Riemann problems can be implemented.
74  */
75 
76  RiemannSolver::RiemannSolver(
78  : m_requiresRotation(false), m_rotStorage (3)
79  {
80  boost::ignore_unused(pSession);
81  }
82 
83  /**
84  * @brief Perform the Riemann solve given the forwards and backwards
85  * spaces.
86  *
87  * This routine calls the virtual function #v_Solve to perform the
88  * Riemann solve. If the flag #m_requiresRotation is set, then the
89  * velocity field is rotated to the normal direction to perform
90  * dimensional splitting, and the resulting fluxes are rotated back to
91  * the Cartesian directions before being returned. For the Rotation to
92  * work, the normal vectors "N" and the location of the vector
93  * components in Fwd "vecLocs"must be set via the SetAuxVec() method.
94  *
95  * @param Fwd Forwards trace space.
96  * @param Bwd Backwards trace space.
97  * @param flux Resultant flux along trace space.
98  */
100  const int nDim,
101  const Array<OneD, const Array<OneD, NekDouble> > &Fwd,
102  const Array<OneD, const Array<OneD, NekDouble> > &Bwd,
104  {
105  if (m_requiresRotation)
106  {
107  ASSERTL1(CheckVectors("N"), "N not defined.");
108  ASSERTL1(CheckAuxVec("vecLocs"), "vecLocs not defined.");
110  m_vectors["N"]();
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] =
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  }
142 
143  /**
144  * @brief Rotate a vector field to trace normal.
145  *
146  * This function performs a rotation of a vector so that the first
147  * component aligns with the trace normal direction.
148  *
149  * The vectors components are stored in inarray. Their locations must
150  * be specified in the "vecLocs" array. vecLocs[0] contains the locations
151  * of the first vectors components, vecLocs[1] those of the second and
152  * so on.
153  *
154  * In 2D, this is accomplished through the transform:
155  *
156  * \f[ (u_x, u_y) = (n_x u_x + n_y u_y, -n_x v_x + n_y v_y) \f]
157  *
158  * In 3D, we generate a (non-unique) transformation using
159  * RiemannSolver::fromToRotation.
160  *
161  */
163  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
164  const Array<OneD, const Array<OneD, NekDouble> > &normals,
165  const Array<OneD, const Array<OneD, NekDouble> > &vecLocs,
166  Array<OneD, Array<OneD, NekDouble> > &outarray)
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  }
244 
245  /**
246  * @brief Rotate a vector field from trace normal.
247  *
248  * This function performs a rotation of the triad of vector components
249  * provided in inarray so that the first component aligns with the
250  * Cartesian components; it performs the inverse operation of
251  * RiemannSolver::rotateToNormal.
252  */
254  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
255  const Array<OneD, const Array<OneD, NekDouble> > &normals,
256  const Array<OneD, const Array<OneD, NekDouble> > &vecLocs,
257  Array<OneD, Array<OneD, NekDouble> > &outarray)
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  }
328 
329  /**
330  * @brief Determine whether a scalar has been defined in #m_scalars.
331  *
332  * @param name Scalar name.
333  */
335  {
336  return m_scalars.find(name) != m_scalars.end();
337  }
338 
339  /**
340  * @brief Determine whether a vector has been defined in #m_vectors.
341  *
342  * @param name Vector name.
343  */
345  {
346  return m_vectors.find(name) != m_vectors.end();
347  }
348 
349  /**
350  * @brief Determine whether a parameter has been defined in #m_params.
351  *
352  * @param name Parameter name.
353  */
355  {
356  return m_params.find(name) != m_params.end();
357  }
358 
359  /**
360  * @brief Determine whether a scalar has been defined in #m_auxScal.
361  *
362  * @param name Scalar name.
363  */
365  {
366  return m_auxScal.find(name) != m_auxScal.end();
367  }
368 
369  /**
370  * @brief Determine whether a vector has been defined in #m_auxVec.
371  *
372  * @param name Vector name.
373  */
375  {
376  return m_auxVec.find(name) != m_auxVec.end();
377  }
378 
379  /**
380  * @brief Generate rotation matrices for 3D expansions.
381  */
383  const Array<OneD, const Array<OneD, NekDouble> > &normals)
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.
394 
395  for (i = 0; i < 9; ++i)
396  {
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  }
414 
415  /**
416  * @brief A function for creating a rotation matrix that rotates a
417  * vector @a from into another vector @a to.
418  *
419  * Authors: Tomas Möller, John Hughes
420  * "Efficiently Building a Matrix to Rotate One Vector to
421  * Another" Journal of Graphics Tools, 4(4):1-4, 1999
422  *
423  * @param from Normalised 3-vector to rotate from.
424  * @param to Normalised 3-vector to rotate to.
425  * @param out Resulting 3x3 rotation matrix (row-major order).
426  */
430  NekDouble *mat)
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  }
512  }
513 }
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 bool CheckScalars(std::string name)
Determine whether a scalar has been defined in m_scalars.
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
STL namespace.
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
#define CROSS(dest, v1, v2)
std::map< std::string, RSScalarFuncType > m_scalars
Map of scalar function types.
void GenerateRotationMatrices(const Array< OneD, const Array< OneD, NekDouble > > &normals)
Generate rotation matrices for 3D expansions.
RiemannSolverFactory & GetRiemannSolverFactory()
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.
std::map< std::string, RSScalarFuncType > m_auxScal
Map of auxiliary scalar function types.
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_rotStorage
Rotation storage.
double NekDouble
SOLVER_UTILS_EXPORT bool CheckAuxScal(std::string name)
Determine whether a scalar has been defined in m_auxScal.
#define EPSILON
Array< OneD, Array< OneD, NekDouble > > m_rotMat
Rotation matrices for each trace quadrature point.
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.
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
bool m_requiresRotation
Indicates whether the Riemann solver requires a rotation to be applied to the velocity fields...
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 DOT(v1, v2)
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
SOLVER_UTILS_EXPORT bool CheckParams(std::string name)
Determine whether a parameter has been defined in m_params.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
std::map< std::string, RSVecFuncType > m_vectors
Map of vector function types.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
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.
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
Provides a generic Factory class.
Definition: NekFactory.hpp:103