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  { \
44  dest[0] = v1[1] * v2[2] - v1[2] * v2[1]; \
45  dest[1] = v1[2] * v2[0] - v1[0] * v2[2]; \
46  dest[2] = v1[0] * v2[1] - v1[1] * v2[0]; \
47  }
48 
49 #define DOT(v1, v2) (v1[0] * v2[0] + v1[1] * v2[1] + v1[2] * v2[2])
50 
51 #define SUB(dest, v1, v2) \
52  { \
53  dest[0] = v1[0] - v2[0]; \
54  dest[1] = v1[1] - v2[1]; \
55  dest[2] = v1[2] - v2[2]; \
56  }
57 
58 using namespace std;
59 
60 namespace Nektar
61 {
62 namespace SolverUtils
63 {
64 /**
65  * Retrieves the singleton instance of the Riemann solver factory.
66  */
68 {
69  static RiemannSolverFactory instance;
70  return instance;
71 }
72 
73 /**
74  * @class RiemannSolver
75  *
76  * @brief The RiemannSolver class provides an abstract interface under
77  * which solvers for various Riemann problems can be implemented.
78  */
79 
80 RiemannSolver::RiemannSolver() : m_requiresRotation(false), m_rotStorage(3)
81 {
82 }
83 
86  : m_requiresRotation(false), m_rotStorage(3)
87 {
88  boost::ignore_unused(pSession);
89 }
90 
91 /**
92  * @brief Perform the Riemann solve given the forwards and backwards
93  * spaces.
94  *
95  * This routine calls the virtual function #v_Solve to perform the
96  * Riemann solve. If the flag #m_requiresRotation is set, then the
97  * velocity field is rotated to the normal direction to perform
98  * dimensional splitting, and the resulting fluxes are rotated back to
99  * the Cartesian directions before being returned. For the Rotation to
100  * work, the normal vectors "N" and the location of the vector
101  * components in Fwd "vecLocs"must be set via the SetAuxVec() method.
102  *
103  * @param Fwd Forwards trace space.
104  * @param Bwd Backwards trace space.
105  * @param flux Resultant flux along trace space.
106  */
107 void RiemannSolver::Solve(const int nDim,
108  const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
109  const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
111 {
112  if (m_requiresRotation)
113  {
114  ASSERTL1(CheckVectors("N"), "N not defined.");
115  ASSERTL1(CheckAuxVec("vecLocs"), "vecLocs not defined.");
117  m_vectors["N"]();
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  {
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]);
139  v_Solve(nDim, m_rotStorage[0], m_rotStorage[1], m_rotStorage[2]);
140  rotateFromNormal(m_rotStorage[2], normals, vecLocs, flux);
141  }
142  else
143  {
144  v_Solve(nDim, Fwd, Bwd, flux);
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, 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  {
217  GenerateRotationMatrices(normals);
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 }
242 
243 /**
244  * @brief Rotate a vector field from trace normal.
245  *
246  * This function performs a rotation of the triad of vector components
247  * provided in inarray so that the first component aligns with the
248  * Cartesian components; it performs the inverse operation of
249  * RiemannSolver::rotateToNormal.
250  */
252  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
253  const Array<OneD, const Array<OneD, NekDouble>> &normals,
254  const Array<OneD, const Array<OneD, NekDouble>> &vecLocs,
255  Array<OneD, Array<OneD, NekDouble>> &outarray)
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 }
319 
320 /**
321  * @brief Determine whether a scalar has been defined in #m_scalars.
322  *
323  * @param name Scalar name.
324  */
326 {
327  return m_scalars.find(name) != m_scalars.end();
328 }
329 
330 /**
331  * @brief Determine whether a vector has been defined in #m_vectors.
332  *
333  * @param name Vector name.
334  */
336 {
337  return m_vectors.find(name) != m_vectors.end();
338 }
339 
340 /**
341  * @brief Determine whether a parameter has been defined in #m_params.
342  *
343  * @param name Parameter name.
344  */
346 {
347  return m_params.find(name) != m_params.end();
348 }
349 
350 /**
351  * @brief Determine whether a scalar has been defined in #m_auxScal.
352  *
353  * @param name Scalar name.
354  */
356 {
357  return m_auxScal.find(name) != m_auxScal.end();
358 }
359 
360 /**
361  * @brief Determine whether a vector has been defined in #m_auxVec.
362  *
363  * @param name Vector name.
364  */
366 {
367  return m_auxVec.find(name) != m_auxVec.end();
368 }
369 
370 /**
371  * @brief Generate rotation matrices for 3D expansions.
372  */
374  const Array<OneD, const Array<OneD, NekDouble>> &normals)
375 {
376  Array<OneD, NekDouble> xdir(3, 0.0);
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.
385 
386  for (i = 0; i < 9; ++i)
387  {
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 }
405 
406 /**
407  * @brief A function for creating a rotation matrix that rotates a
408  * vector @a from into another vector @a to.
409  *
410  * Authors: Tomas Möller, John Hughes
411  * "Efficiently Building a Matrix to Rotate One Vector to
412  * Another" Journal of Graphics Tools, 4(4):1-4, 1999
413  *
414  * @param from Normalised 3-vector to rotate from.
415  * @param to Normalised 3-vector to rotate to.
416  * @param out Resulting 3x3 rotation matrix (row-major order).
417  */
420  NekDouble *mat)
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 }
507 
508 /**
509  * @brief Calculate the flux jacobian of Fwd and Bwd
510  *
511  * @param Fwd Forwards trace space.
512  * @param Bwd Backwards trace space.
513  * @param flux Resultant flux along trace space.
514  */
516  const int nDim, const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
517  const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
519 {
520  int nPts = Fwd[0].size();
521 
522  if (m_requiresRotation)
523  {
524  ASSERTL1(CheckVectors("N"), "N not defined.");
525  ASSERTL1(CheckAuxVec("vecLocs"), "vecLocs not defined.");
527  m_vectors["N"]();
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 }
545 
547  const int nDim, const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
548  const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
549  const Array<OneD, const Array<OneD, NekDouble>> &normals,
551 {
552  boost::ignore_unused(nDim, Fwd, Bwd, normals, FJac, BJac);
553  NEKERROR(ErrorUtil::efatal, "v_CalcFluxJacobian not specified.");
554 }
555 
556 } // namespace SolverUtils
557 } // namespace Nektar
#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:249
#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 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.
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.
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.
SOLVER_UTILS_EXPORT void GenerateRotationMatrices(const Array< OneD, const Array< OneD, NekDouble >> &normals)
Generate rotation matrices for 3D expansions.
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 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.
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.
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)
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.
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.
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 bool CheckAuxScal(std::string name)
Determine whether a scalar has been defined in m_auxScal.
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
SOLVER_UTILS_EXPORT RiemannSolver()
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:2
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:77
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:209
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:574
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:598
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:692
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255