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
37
38#define EPSILON 0.000001
39
40#define CROSS(dest, v1, v2) \
41 { \
42 dest[0] = v1[1] * v2[2] - v1[2] * v2[1]; \
43 dest[1] = v1[2] * v2[0] - v1[0] * v2[2]; \
44 dest[2] = v1[0] * v2[1] - v1[1] * v2[0]; \
45 }
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 { \
51 dest[0] = v1[0] - v2[0]; \
52 dest[1] = v1[1] - v2[1]; \
53 dest[2] = v1[2] - v2[2]; \
54 }
55
56using namespace std;
57
58namespace Nektar::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
76RiemannSolver::RiemannSolver() : m_requiresRotation(false), m_rotStorage(3)
77{
78}
79
81 [[maybe_unused]] const LibUtilities::SessionReaderSharedPtr &pSession)
82 : m_requiresRotation(false), m_rotStorage(3)
83{
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 */
102void RiemannSolver::Solve(const int nDim,
103 const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
104 const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
106{
108 {
109 ASSERTL1(CheckVectors("N"), "N not defined.");
110 ASSERTL1(CheckAuxVec("vecLocs"), "vecLocs not defined.");
112 m_vectors["N"]();
114 m_auxVec["vecLocs"]();
115
116 int nFields = Fwd.size();
117 int nPts = Fwd[0].size();
118
119 if (m_rotStorage[0].size() != nFields ||
120 m_rotStorage[0][0].size() != nPts)
121 {
122 for (int i = 0; i < 3; ++i)
123 {
125 for (int j = 0; j < nFields; ++j)
126 {
128 }
129 }
130 }
131
132 rotateToNormal(Fwd, normals, vecLocs, m_rotStorage[0]);
133 rotateToNormal(Bwd, normals, vecLocs, m_rotStorage[1]);
134
136 rotateFromNormal(m_rotStorage[2], normals, vecLocs, flux);
137
138 // Attempt to subtract (\vec{U}\vec{vg})\dot n for ALE
139 if (m_ALESolver)
140 {
141 auto vgt = m_vectors["vgt"]();
142 auto N = m_vectors["N"]();
143 for (int i = 0; i < flux.size(); ++i)
144 {
145 for (int j = 0; j < flux[i].size(); ++j)
146 {
147 NekDouble tmp = 0;
148 for (int k = 0; k < nDim; ++k)
149 {
150 tmp += N[k][j] * vgt[k][j];
151 }
152 flux[i][j] -= 0.5 * (Fwd[i][j] + Bwd[i][j]) * tmp;
153 }
154 }
155 }
156 }
157 else
158 {
159 v_Solve(nDim, Fwd, Bwd, flux);
160 }
161}
162
163/**
164 * @brief Rotate a vector field to trace normal.
165 *
166 * This function performs a rotation of a vector so that the first
167 * component aligns with the trace normal direction.
168 *
169 * The vectors components are stored in inarray. Their locations must
170 * be specified in the "vecLocs" array. vecLocs[0] contains the locations
171 * of the first vectors components, vecLocs[1] those of the second and
172 * so on.
173 *
174 * In 2D, this is accomplished through the transform:
175 *
176 * \f[ (u_x, u_y) = (n_x u_x + n_y u_y, -n_x v_x + n_y v_y) \f]
177 *
178 * In 3D, we generate a (non-unique) transformation using
179 * RiemannSolver::fromToRotation.
180 *
181 */
183 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
184 const Array<OneD, const Array<OneD, NekDouble>> &normals,
185 const Array<OneD, const Array<OneD, NekDouble>> &vecLocs,
187{
188 for (int i = 0; i < inarray.size(); ++i)
189 {
190 Vmath::Vcopy(inarray[i].size(), inarray[i], 1, outarray[i], 1);
191 }
192
193 for (int i = 0; i < vecLocs.size(); i++)
194 {
195 ASSERTL1(vecLocs[i].size() == normals.size(),
196 "vecLocs[i] element count mismatch");
197
198 switch (normals.size())
199 {
200 case 1:
201 { // do nothing
202 const int nq = inarray[0].size();
203 const int vx = (int)vecLocs[i][0];
204 Vmath::Vmul(nq, inarray[vx], 1, normals[0], 1, outarray[vx], 1);
205 break;
206 }
207 case 2:
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
213 Vmath::Vmul(nq, inarray[vx], 1, normals[0], 1, outarray[vx], 1);
214 Vmath::Vvtvp(nq, inarray[vy], 1, normals[1], 1, outarray[vx], 1,
215 outarray[vx], 1);
216 Vmath::Vmul(nq, inarray[vx], 1, normals[1], 1, outarray[vy], 1);
217 Vmath::Vvtvm(nq, inarray[vy], 1, normals[0], 1, outarray[vy], 1,
218 outarray[vy], 1);
219 break;
220 }
221
222 case 3:
223 {
224 const int nq = inarray[0].size();
225 const int vx = (int)vecLocs[i][0];
226 const int vy = (int)vecLocs[i][1];
227 const int vz = (int)vecLocs[i][2];
228
229 // Generate matrices if they don't already exist.
230 if (m_rotMat.size() == 0)
231 {
233 }
234
235 // Apply rotation matrices.
236 Vmath::Vvtvvtp(nq, inarray[vx], 1, m_rotMat[0], 1, inarray[vy],
237 1, m_rotMat[1], 1, outarray[vx], 1);
238 Vmath::Vvtvp(nq, inarray[vz], 1, m_rotMat[2], 1, outarray[vx],
239 1, outarray[vx], 1);
240 Vmath::Vvtvvtp(nq, inarray[vx], 1, m_rotMat[3], 1, inarray[vy],
241 1, m_rotMat[4], 1, outarray[vy], 1);
242 Vmath::Vvtvp(nq, inarray[vz], 1, m_rotMat[5], 1, outarray[vy],
243 1, outarray[vy], 1);
244 Vmath::Vvtvvtp(nq, inarray[vx], 1, m_rotMat[6], 1, inarray[vy],
245 1, m_rotMat[7], 1, outarray[vz], 1);
246 Vmath::Vvtvp(nq, inarray[vz], 1, m_rotMat[8], 1, outarray[vz],
247 1, outarray[vz], 1);
248 break;
249 }
250
251 default:
252 ASSERTL1(false, "Invalid space dimension.");
253 break;
254 }
255 }
256}
257
258/**
259 * @brief Rotate a vector field from trace normal.
260 *
261 * This function performs a rotation of the triad of vector components
262 * provided in inarray so that the first component aligns with the
263 * Cartesian components; it performs the inverse operation of
264 * RiemannSolver::rotateToNormal.
265 */
267 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
268 const Array<OneD, const Array<OneD, NekDouble>> &normals,
269 const Array<OneD, const Array<OneD, NekDouble>> &vecLocs,
271{
272 for (int i = 0; i < inarray.size(); ++i)
273 {
274 Vmath::Vcopy(inarray[i].size(), inarray[i], 1, outarray[i], 1);
275 }
276
277 for (int i = 0; i < vecLocs.size(); i++)
278 {
279 ASSERTL1(vecLocs[i].size() == normals.size(),
280 "vecLocs[i] element count mismatch");
281
282 switch (normals.size())
283 {
284 case 1:
285 { // do nothing
286 const int nq = normals[0].size();
287 const int vx = (int)vecLocs[i][0];
288 Vmath::Vmul(nq, inarray[vx], 1, normals[0], 1, outarray[vx], 1);
289 break;
290 }
291 case 2:
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
297 Vmath::Vmul(nq, inarray[vy], 1, normals[1], 1, outarray[vx], 1);
298 Vmath::Vvtvm(nq, inarray[vx], 1, normals[0], 1, outarray[vx], 1,
299 outarray[vx], 1);
300 Vmath::Vmul(nq, inarray[vx], 1, normals[1], 1, outarray[vy], 1);
301 Vmath::Vvtvp(nq, inarray[vy], 1, normals[0], 1, outarray[vy], 1,
302 outarray[vy], 1);
303 break;
304 }
305
306 case 3:
307 {
308 const int nq = normals[0].size();
309 const int vx = (int)vecLocs[i][0];
310 const int vy = (int)vecLocs[i][1];
311 const int vz = (int)vecLocs[i][2];
312
313 Vmath::Vvtvvtp(nq, inarray[vx], 1, m_rotMat[0], 1, inarray[vy],
314 1, m_rotMat[3], 1, outarray[vx], 1);
315 Vmath::Vvtvp(nq, inarray[vz], 1, m_rotMat[6], 1, outarray[vx],
316 1, outarray[vx], 1);
317 Vmath::Vvtvvtp(nq, inarray[vx], 1, m_rotMat[1], 1, inarray[vy],
318 1, m_rotMat[4], 1, outarray[vy], 1);
319 Vmath::Vvtvp(nq, inarray[vz], 1, m_rotMat[7], 1, outarray[vy],
320 1, outarray[vy], 1);
321 Vmath::Vvtvvtp(nq, inarray[vx], 1, m_rotMat[2], 1, inarray[vy],
322 1, m_rotMat[5], 1, outarray[vz], 1);
323 Vmath::Vvtvp(nq, inarray[vz], 1, m_rotMat[8], 1, outarray[vz],
324 1, outarray[vz], 1);
325 break;
326 }
327
328 default:
329 ASSERTL1(false, "Invalid space dimension.");
330 break;
331 }
332 }
333}
334
335/**
336 * @brief Determine whether a scalar has been defined in #m_scalars.
337 *
338 * @param name Scalar name.
339 */
341{
342 return m_scalars.find(name) != m_scalars.end();
343}
344
345/**
346 * @brief Determine whether a vector has been defined in #m_vectors.
347 *
348 * @param name Vector name.
349 */
351{
352 return m_vectors.find(name) != m_vectors.end();
353}
354
355/**
356 * @brief Determine whether a parameter has been defined in #m_params.
357 *
358 * @param name Parameter name.
359 */
361{
362 return m_params.find(name) != m_params.end();
363}
364
365/**
366 * @brief Determine whether a scalar has been defined in #m_auxScal.
367 *
368 * @param name Scalar name.
369 */
371{
372 return m_auxScal.find(name) != m_auxScal.end();
373}
374
375/**
376 * @brief Determine whether a vector has been defined in #m_auxVec.
377 *
378 * @param name Vector name.
379 */
381{
382 return m_auxVec.find(name) != m_auxVec.end();
383}
384
385/**
386 * @brief Generate rotation matrices for 3D expansions.
387 */
389 const Array<OneD, const Array<OneD, NekDouble>> &normals)
390{
391 Array<OneD, NekDouble> xdir(3, 0.0);
393 NekDouble tmp[9];
394 const int nq = normals[0].size();
395 int i, j;
396 xdir[0] = 1.0;
397
398 // Allocate storage for rotation matrices.
400
401 for (i = 0; i < 9; ++i)
402 {
404 }
405 for (i = 0; i < normals[0].size(); ++i)
406 {
407 // Generate matrix which takes us from (1,0,0) vector to trace
408 // normal.
409 tn[0] = normals[0][i];
410 tn[1] = normals[1][i];
411 tn[2] = normals[2][i];
412 FromToRotation(tn, xdir, tmp);
413
414 for (j = 0; j < 9; ++j)
415 {
416 m_rotMat[j][i] = tmp[j];
417 }
418 }
419}
420
421/**
422 * @brief A function for creating a rotation matrix that rotates a
423 * vector @a from into another vector @a to.
424 *
425 * Authors: Tomas Möller, John Hughes
426 * "Efficiently Building a Matrix to Rotate One Vector to
427 * Another" Journal of Graphics Tools, 4(4):1-4, 1999
428 *
429 * @param from Normalised 3-vector to rotate from.
430 * @param to Normalised 3-vector to rotate to.
431 * @param out Resulting 3x3 rotation matrix (row-major order).
432 */
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;
459 x[1] = x[2] = 0.0;
460 }
461 else
462 {
463 x[2] = 1.0;
464 x[0] = x[1] = 0.0;
465 }
466 }
467 else
468 {
469 if (x[1] < x[2])
470 {
471 x[1] = 1.0;
472 x[0] = x[2] = 0.0;
473 }
474 else
475 {
476 x[2] = 1.0;
477 x[0] = x[1] = 0.0;
478 }
479 }
480
481 u[0] = x[0] - from[0];
482 u[1] = x[1] - from[1];
483 u[2] = x[2] - from[2];
484 v[0] = x[0] - to[0];
485 v[1] = x[1] - to[1];
486 v[2] = x[2] - to[2];
487
488 c1 = 2.0 / DOT(u, u);
489 c2 = 2.0 / DOT(v, v);
490 c3 = c1 * c2 * DOT(u, v);
491
492 for (i = 0; i < 3; i++)
493 {
494 for (j = 0; j < 3; j++)
495 {
496 mat[3 * i + j] =
497 -c1 * u[i] * u[j] - c2 * v[i] * v[j] + c3 * v[i] * u[j];
498 }
499 mat[i + 3 * i] += 1.0;
500 }
501 }
502 else
503 {
504 NekDouble hvx, hvz, hvxy, hvxz, hvyz;
505 h = 1.0 / (1.0 + e);
506 hvx = h * v[0];
507 hvz = h * v[2];
508 hvxy = hvx * v[1];
509 hvxz = hvx * v[2];
510 hvyz = hvz * v[1];
511 mat[0] = e + hvx * v[0];
512 mat[1] = hvxy - v[2];
513 mat[2] = hvxz + v[1];
514 mat[3] = hvxy + v[2];
515 mat[4] = e + h * v[1] * v[1];
516 mat[5] = hvyz - v[0];
517 mat[6] = hvxz - v[1];
518 mat[7] = hvyz + v[0];
519 mat[8] = e + hvz * v[2];
520 }
521}
522
523/**
524 * @brief Calculate the flux jacobian of Fwd and Bwd
525 *
526 * @param Fwd Forwards trace space.
527 * @param Bwd Backwards trace space.
528 * @param flux Resultant flux along trace space.
529 */
531 const int nDim, const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
532 const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
534{
535 int nPts = Fwd[0].size();
536
538 {
539 ASSERTL1(CheckVectors("N"), "N not defined.");
540 ASSERTL1(CheckAuxVec("vecLocs"), "vecLocs not defined.");
542 m_vectors["N"]();
544 m_auxVec["vecLocs"]();
545
546 v_CalcFluxJacobian(nDim, Fwd, Bwd, normals, FJac, BJac);
547 }
548 else
549 {
551 for (int i = 0; i < nDim; i++)
552 {
553 normals[i] = Array<OneD, NekDouble>(nPts, 0.0);
554 }
555 Vmath::Fill(nPts, 1.0, normals[0], 1);
556
557 v_CalcFluxJacobian(nDim, Fwd, Bwd, normals, FJac, BJac);
558 }
559}
560
562 [[maybe_unused]] const int nDim,
563 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
564 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
565 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>> &normals,
566 [[maybe_unused]] DNekBlkMatSharedPtr &FJac,
567 [[maybe_unused]] DNekBlkMatSharedPtr &BJac)
568{
569 NEKERROR(ErrorUtil::efatal, "v_CalcFluxJacobian not specified.");
570}
571
572} // namespace Nektar::SolverUtils
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
#define EPSILON
#define CROSS(dest, v1, v2)
#define DOT(v1, v2)
Provides a generic Factory class.
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.
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
SOLVER_UTILS_EXPORT RiemannSolver()
bool m_ALESolver
Flag if using the ALE formulation.
std::map< std::string, RSVecFuncType > m_auxVec
Map of auxiliary vector function types.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
RiemannSolverFactory & GetRiemannSolverFactory()
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.hpp:72
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.hpp:366
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.hpp:381
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.hpp:54
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.hpp:439
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
STL namespace.