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()
77  : m_requiresRotation(false), m_rotStorage(3) {}
78 
81  : m_requiresRotation(false), m_rotStorage(3)
82  {
83  boost::ignore_unused(pSession);
84  }
85 
86  /**
87  * @brief Perform the Riemann solve given the forwards and backwards
88  * spaces.
89  *
90  * This routine calls the virtual function #v_Solve to perform the
91  * Riemann solve. If the flag #m_requiresRotation is set, then the
92  * velocity field is rotated to the normal direction to perform
93  * dimensional splitting, and the resulting fluxes are rotated back to
94  * the Cartesian directions before being returned. For the Rotation to
95  * work, the normal vectors "N" and the location of the vector
96  * components in Fwd "vecLocs"must be set via the SetAuxVec() method.
97  *
98  * @param Fwd Forwards trace space.
99  * @param Bwd Backwards trace space.
100  * @param flux Resultant flux along trace space.
101  */
103  const int nDim,
104  const Array<OneD, const Array<OneD, NekDouble> > &Fwd,
105  const Array<OneD, const Array<OneD, NekDouble> > &Bwd,
107  {
108  if (m_requiresRotation)
109  {
110  ASSERTL1(CheckVectors("N"), "N not defined.");
111  ASSERTL1(CheckAuxVec("vecLocs"), "vecLocs not defined.");
113  m_vectors["N"]();
115  m_auxVec["vecLocs"]();
116 
117  int nFields = Fwd .size();
118  int nPts = Fwd[0].size();
119 
120  if (m_rotStorage[0].size() != nFields ||
121  m_rotStorage[0][0].size() != nPts)
122  {
123  for (int i = 0; i < 3; ++i)
124  {
125  m_rotStorage[i] =
127  for (int j = 0; j < nFields; ++j)
128  {
129  m_rotStorage[i][j] = Array<OneD, NekDouble>(nPts);
130  }
131  }
132  }
133 
134  rotateToNormal (Fwd, normals, vecLocs, m_rotStorage[0]);
135  rotateToNormal (Bwd, normals, vecLocs, m_rotStorage[1]);
136  v_Solve (nDim, m_rotStorage[0], m_rotStorage[1],
137  m_rotStorage[2]);
138  rotateFromNormal(m_rotStorage[2], normals, vecLocs, flux);
139  }
140  else
141  {
142  v_Solve(nDim, Fwd, Bwd, flux);
143  }
144  }
145 
146 
147 
148  /**
149  * @brief Rotate a vector field to trace normal.
150  *
151  * This function performs a rotation of a vector so that the first
152  * component aligns with the trace normal direction.
153  *
154  * The vectors components are stored in inarray. Their locations must
155  * be specified in the "vecLocs" array. vecLocs[0] contains the locations
156  * of the first vectors components, vecLocs[1] those of the second and
157  * so on.
158  *
159  * In 2D, this is accomplished through the transform:
160  *
161  * \f[ (u_x, u_y) = (n_x u_x + n_y u_y, -n_x v_x + n_y v_y) \f]
162  *
163  * In 3D, we generate a (non-unique) transformation using
164  * RiemannSolver::fromToRotation.
165  *
166  */
168  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
169  const Array<OneD, const Array<OneD, NekDouble> > &normals,
170  const Array<OneD, const Array<OneD, NekDouble> > &vecLocs,
171  Array<OneD, Array<OneD, NekDouble> > &outarray)
172  {
173  for (int i = 0; i < inarray.size(); ++i)
174  {
175  Vmath::Vcopy(inarray[i].size(), inarray[i], 1,
176  outarray[i], 1);
177  }
178 
179  for (int i = 0; i < vecLocs.size(); i++)
180  {
181  ASSERTL1(vecLocs[i].size() == normals.size(),
182  "vecLocs[i] element count mismatch");
183 
184  switch (normals.size())
185  {
186  case 1:
187  { // do nothing
188  const int nq = inarray[0].size();
189  const int vx = (int)vecLocs[i][0];
190  Vmath::Vmul (nq, inarray [vx], 1, normals [0], 1,
191  outarray[vx], 1);
192  break;
193  }
194  case 2:
195  {
196  const int nq = inarray[0].size();
197  const int vx = (int)vecLocs[i][0];
198  const int vy = (int)vecLocs[i][1];
199 
200  Vmath::Vmul (nq, inarray [vx], 1, normals [0], 1,
201  outarray[vx], 1);
202  Vmath::Vvtvp(nq, inarray [vy], 1, normals [1], 1,
203  outarray[vx], 1, outarray[vx], 1);
204  Vmath::Vmul (nq, inarray [vx], 1, normals [1], 1,
205  outarray[vy], 1);
206  Vmath::Vvtvm(nq, inarray [vy], 1, normals [0], 1,
207  outarray[vy], 1, outarray[vy], 1);
208  break;
209  }
210 
211  case 3:
212  {
213  const int nq = inarray[0].size();
214  const int vx = (int)vecLocs[i][0];
215  const int vy = (int)vecLocs[i][1];
216  const int vz = (int)vecLocs[i][2];
217 
218  // Generate matrices if they don't already exist.
219  if (m_rotMat.size() == 0)
220  {
221  GenerateRotationMatrices(normals);
222  }
223 
224  // Apply rotation matrices.
225  Vmath::Vvtvvtp(nq, inarray [vx], 1, m_rotMat[0], 1,
226  inarray [vy], 1, m_rotMat[1], 1,
227  outarray[vx], 1);
228  Vmath::Vvtvp (nq, inarray [vz], 1, m_rotMat[2], 1,
229  outarray[vx], 1, outarray[vx], 1);
230  Vmath::Vvtvvtp(nq, inarray [vx], 1, m_rotMat[3], 1,
231  inarray [vy], 1, m_rotMat[4], 1,
232  outarray[vy], 1);
233  Vmath::Vvtvp (nq, inarray [vz], 1, m_rotMat[5], 1,
234  outarray[vy], 1, outarray[vy], 1);
235  Vmath::Vvtvvtp(nq, inarray [vx], 1, m_rotMat[6], 1,
236  inarray [vy], 1, m_rotMat[7], 1,
237  outarray[vz], 1);
238  Vmath::Vvtvp (nq, inarray [vz], 1, m_rotMat[8], 1,
239  outarray[vz], 1, outarray[vz], 1);
240  break;
241  }
242 
243  default:
244  ASSERTL1(false, "Invalid space dimension.");
245  break;
246  }
247  }
248  }
249 
250  /**
251  * @brief Rotate a vector field from trace normal.
252  *
253  * This function performs a rotation of the triad of vector components
254  * provided in inarray so that the first component aligns with the
255  * Cartesian components; it performs the inverse operation of
256  * RiemannSolver::rotateToNormal.
257  */
259  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
260  const Array<OneD, const Array<OneD, NekDouble> > &normals,
261  const Array<OneD, const Array<OneD, NekDouble> > &vecLocs,
262  Array<OneD, Array<OneD, NekDouble> > &outarray)
263  {
264  for (int i = 0; i < inarray.size(); ++i)
265  {
266  Vmath::Vcopy(inarray[i].size(), inarray[i], 1,
267  outarray[i], 1);
268  }
269 
270  for (int i = 0; i < vecLocs.size(); i++)
271  {
272  ASSERTL1(vecLocs[i].size() == normals.size(),
273  "vecLocs[i] element count mismatch");
274 
275  switch (normals.size())
276  {
277  case 1:
278  { // do nothing
279  const int nq = normals[0].size();
280  const int vx = (int)vecLocs[i][0];
281  Vmath::Vmul (nq, inarray [vx], 1, normals [0], 1,
282  outarray[vx], 1);
283  break;
284  }
285  case 2:
286  {
287  const int nq = normals[0].size();
288  const int vx = (int)vecLocs[i][0];
289  const int vy = (int)vecLocs[i][1];
290 
291  Vmath::Vmul (nq, inarray [vy], 1, normals [1], 1,
292  outarray[vx], 1);
293  Vmath::Vvtvm(nq, inarray [vx], 1, normals [0], 1,
294  outarray[vx], 1, outarray[vx], 1);
295  Vmath::Vmul (nq, inarray [vx], 1, normals [1], 1,
296  outarray[vy], 1);
297  Vmath::Vvtvp(nq, inarray [vy], 1, normals [0], 1,
298  outarray[vy], 1, outarray[vy], 1);
299  break;
300  }
301 
302  case 3:
303  {
304  const int nq = normals[0].size();
305  const int vx = (int)vecLocs[i][0];
306  const int vy = (int)vecLocs[i][1];
307  const int vz = (int)vecLocs[i][2];
308 
309  Vmath::Vvtvvtp(nq, inarray [vx], 1, m_rotMat[0], 1,
310  inarray [vy], 1, m_rotMat[3], 1,
311  outarray[vx], 1);
312  Vmath::Vvtvp (nq, inarray [vz], 1, m_rotMat[6], 1,
313  outarray[vx], 1, outarray[vx], 1);
314  Vmath::Vvtvvtp(nq, inarray [vx], 1, m_rotMat[1], 1,
315  inarray [vy], 1, m_rotMat[4], 1,
316  outarray[vy], 1);
317  Vmath::Vvtvp (nq, inarray [vz], 1, m_rotMat[7], 1,
318  outarray[vy], 1, outarray[vy], 1);
319  Vmath::Vvtvvtp(nq, inarray [vx], 1, m_rotMat[2], 1,
320  inarray [vy], 1, m_rotMat[5], 1,
321  outarray[vz], 1);
322  Vmath::Vvtvp (nq, inarray [vz], 1, m_rotMat[8], 1,
323  outarray[vz], 1, outarray[vz], 1);
324  break;
325  }
326 
327  default:
328  ASSERTL1(false, "Invalid space dimension.");
329  break;
330  }
331  }
332  }
333 
334  /**
335  * @brief Determine whether a scalar has been defined in #m_scalars.
336  *
337  * @param name Scalar name.
338  */
340  {
341  return m_scalars.find(name) != m_scalars.end();
342  }
343 
344  /**
345  * @brief Determine whether a vector has been defined in #m_vectors.
346  *
347  * @param name Vector name.
348  */
350  {
351  return m_vectors.find(name) != m_vectors.end();
352  }
353 
354  /**
355  * @brief Determine whether a parameter has been defined in #m_params.
356  *
357  * @param name Parameter name.
358  */
360  {
361  return m_params.find(name) != m_params.end();
362  }
363 
364  /**
365  * @brief Determine whether a scalar has been defined in #m_auxScal.
366  *
367  * @param name Scalar name.
368  */
370  {
371  return m_auxScal.find(name) != m_auxScal.end();
372  }
373 
374  /**
375  * @brief Determine whether a vector has been defined in #m_auxVec.
376  *
377  * @param name Vector name.
378  */
380  {
381  return m_auxVec.find(name) != m_auxVec.end();
382  }
383 
384  /**
385  * @brief Generate rotation matrices for 3D expansions.
386  */
388  const Array<OneD, const Array<OneD, NekDouble> > &normals)
389  {
390  Array<OneD, NekDouble> xdir(3,0.0);
391  Array<OneD, NekDouble> tn (3);
392  NekDouble tmp[9];
393  const int nq = normals[0].size();
394  int i, j;
395  xdir[0] = 1.0;
396 
397  // Allocate storage for rotation matrices.
399 
400  for (i = 0; i < 9; ++i)
401  {
403  }
404  for (i = 0; i < normals[0].size(); ++i)
405  {
406  // Generate matrix which takes us from (1,0,0) vector to trace
407  // normal.
408  tn[0] = normals[0][i];
409  tn[1] = normals[1][i];
410  tn[2] = normals[2][i];
411  FromToRotation(tn, xdir, tmp);
412 
413  for (j = 0; j < 9; ++j)
414  {
415  m_rotMat[j][i] = tmp[j];
416  }
417  }
418  }
419 
420  /**
421  * @brief A function for creating a rotation matrix that rotates a
422  * vector @a from into another vector @a to.
423  *
424  * Authors: Tomas Möller, John Hughes
425  * "Efficiently Building a Matrix to Rotate One Vector to
426  * Another" Journal of Graphics Tools, 4(4):1-4, 1999
427  *
428  * @param from Normalised 3-vector to rotate from.
429  * @param to Normalised 3-vector to rotate to.
430  * @param out Resulting 3x3 rotation matrix (row-major order).
431  */
435  NekDouble *mat)
436  {
437  NekDouble v[3];
438  NekDouble e, h, f;
439 
440  CROSS(v, from, to);
441  e = DOT(from, to);
442  f = (e < 0)? -e:e;
443  if (f > 1.0 - EPSILON)
444  {
445  NekDouble u[3], v[3];
446  NekDouble x[3];
447  NekDouble c1, c2, c3;
448  int i, j;
449 
450  x[0] = (from[0] > 0.0)? from[0] : -from[0];
451  x[1] = (from[1] > 0.0)? from[1] : -from[1];
452  x[2] = (from[2] > 0.0)? from[2] : -from[2];
453 
454  if (x[0] < x[1])
455  {
456  if (x[0] < x[2])
457  {
458  x[0] = 1.0; x[1] = x[2] = 0.0;
459  }
460  else
461  {
462  x[2] = 1.0; x[0] = x[1] = 0.0;
463  }
464  }
465  else
466  {
467  if (x[1] < x[2])
468  {
469  x[1] = 1.0; x[0] = x[2] = 0.0;
470  }
471  else
472  {
473  x[2] = 1.0; x[0] = x[1] = 0.0;
474  }
475  }
476 
477  u[0] = x[0] - from[0];
478  u[1] = x[1] - from[1];
479  u[2] = x[2] - from[2];
480  v[0] = x[0] - to [0];
481  v[1] = x[1] - to [1];
482  v[2] = x[2] - to [2];
483 
484  c1 = 2.0 / DOT(u, u);
485  c2 = 2.0 / DOT(v, v);
486  c3 = c1 * c2 * DOT(u, v);
487 
488  for (i = 0; i < 3; i++) {
489  for (j = 0; j < 3; j++) {
490  mat[3*i+j] = - c1 * u[i] * u[j]
491  - c2 * v[i] * v[j]
492  + c3 * v[i] * u[j];
493  }
494  mat[i+3*i] += 1.0;
495  }
496  }
497  else
498  {
499  NekDouble hvx, hvz, hvxy, hvxz, hvyz;
500  h = 1.0/(1.0 + e);
501  hvx = h * v[0];
502  hvz = h * v[2];
503  hvxy = hvx * v[1];
504  hvxz = hvx * v[2];
505  hvyz = hvz * v[1];
506  mat[0] = e + hvx * v[0];
507  mat[1] = hvxy - v[2];
508  mat[2] = hvxz + v[1];
509  mat[3] = hvxy + v[2];
510  mat[4] = e + h * v[1] * v[1];
511  mat[5] = hvyz - v[0];
512  mat[6] = hvxz - v[1];
513  mat[7] = hvyz + v[0];
514  mat[8] = e + hvz * v[2];
515  }
516  }
517 
518  /**
519  * @brief Calculate the flux jacobian of Fwd and Bwd
520  *
521  * @param Fwd Forwards trace space.
522  * @param Bwd Backwards trace space.
523  * @param flux Resultant flux along trace space.
524  */
526  const int nDim,
527  const Array<OneD, const Array<OneD, NekDouble> > &Fwd,
528  const Array<OneD, const Array<OneD, NekDouble> > &Bwd,
529  DNekBlkMatSharedPtr &FJac,
530  DNekBlkMatSharedPtr &BJac)
531  {
532  int nPts = Fwd[0].size();
533 
534  if (m_requiresRotation)
535  {
536  ASSERTL1(CheckVectors("N"), "N not defined.");
537  ASSERTL1(CheckAuxVec("vecLocs"), "vecLocs not defined.");
539  m_vectors["N"]();
541  m_auxVec["vecLocs"]();
542 
543  v_CalcFluxJacobian(nDim, Fwd, Bwd,normals, FJac, BJac);
544  }
545  else
546  {
547  Array<OneD, Array<OneD, NekDouble> > normals(nDim);
548  for(int i=0;i< nDim;i++)
549  {
550  normals[i] = Array<OneD, NekDouble> (nPts,0.0);
551  }
552  Vmath::Fill(nPts, 1.0, normals[0],1);
553 
554  v_CalcFluxJacobian(nDim, Fwd, Bwd,normals, FJac, BJac);
555  }
556  }
557 
558 
560  const int nDim,
561  const Array<OneD, const Array<OneD, NekDouble> > &Fwd,
562  const Array<OneD, const Array<OneD, NekDouble> > &Bwd,
563  const Array<OneD, const Array<OneD, NekDouble> > &normals,
564  DNekBlkMatSharedPtr &FJac,
565  DNekBlkMatSharedPtr &BJac)
566  {
567  boost::ignore_unused(nDim,Fwd,Bwd,normals,FJac,BJac);
568  NEKERROR(ErrorUtil::efatal, "v_CalcFluxJacobian not specified.");
569  }
570 
571  }
572 }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
#define EPSILON
#define CROSS(dest, v1, v2)
#define DOT(v1, v2)
Provides a generic Factory class.
Definition: NekFactory.hpp:105
SOLVER_UTILS_EXPORT bool CheckAuxVec(std::string name)
Determine whether a vector has been defined in m_auxVec.
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.
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.
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)
SOLVER_UTILS_EXPORT bool CheckParams(std::string name)
Determine whether a parameter has been defined in m_params.
std::map< std::string, RSScalarFuncType > m_auxScal
Map of auxiliary scalar 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.
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.
Array< OneD, Array< OneD, NekDouble > > m_rotMat
Rotation matrices for each trace quadrature point.
SOLVER_UTILS_EXPORT bool CheckScalars(std::string name)
Determine whether a scalar has been defined in m_scalars.
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.
std::map< std::string, RSVecFuncType > m_vectors
Map of vector function types.
std::map< std::string, RSScalarFuncType > m_scalars
Map of scalar function types.
SOLVER_UTILS_EXPORT bool CheckVectors(std::string name)
Determine whether a vector has been defined in m_vectors.
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 GenerateRotationMatrices(const Array< OneD, const Array< OneD, NekDouble > > &normals)
Generate rotation matrices for 3D expansions.
SOLVER_UTILS_EXPORT bool CheckAuxScal(std::string name)
Determine whether a scalar has been defined in m_auxScal.
SOLVER_UTILS_EXPORT RiemannSolver()
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
std::map< std::string, RSVecFuncType > m_auxVec
Map of auxiliary vector function types.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
RiemannSolverFactory & GetRiemannSolverFactory()
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:71
double NekDouble
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:192
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:513
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:541
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45
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:625
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199