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]);
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,
167{
168 for (int i = 0; i < inarray.size(); ++i)
169 {
170 Vmath::Vcopy(inarray[i].size(), inarray[i], 1, outarray[i], 1);
171 }
172
173 for (int i = 0; i < vecLocs.size(); i++)
174 {
175 ASSERTL1(vecLocs[i].size() == normals.size(),
176 "vecLocs[i] element count mismatch");
177
178 switch (normals.size())
179 {
180 case 1:
181 { // do nothing
182 const int nq = inarray[0].size();
183 const int vx = (int)vecLocs[i][0];
184 Vmath::Vmul(nq, inarray[vx], 1, normals[0], 1, outarray[vx], 1);
185 break;
186 }
187 case 2:
188 {
189 const int nq = inarray[0].size();
190 const int vx = (int)vecLocs[i][0];
191 const int vy = (int)vecLocs[i][1];
192
193 Vmath::Vmul(nq, inarray[vx], 1, normals[0], 1, outarray[vx], 1);
194 Vmath::Vvtvp(nq, inarray[vy], 1, normals[1], 1, outarray[vx], 1,
195 outarray[vx], 1);
196 Vmath::Vmul(nq, inarray[vx], 1, normals[1], 1, outarray[vy], 1);
197 Vmath::Vvtvm(nq, inarray[vy], 1, normals[0], 1, outarray[vy], 1,
198 outarray[vy], 1);
199 break;
200 }
201
202 case 3:
203 {
204 const int nq = inarray[0].size();
205 const int vx = (int)vecLocs[i][0];
206 const int vy = (int)vecLocs[i][1];
207 const int vz = (int)vecLocs[i][2];
208
209 // Generate matrices if they don't already exist.
210 if (m_rotMat.size() == 0)
211 {
213 }
214
215 // Apply rotation matrices.
216 Vmath::Vvtvvtp(nq, inarray[vx], 1, m_rotMat[0], 1, inarray[vy],
217 1, m_rotMat[1], 1, outarray[vx], 1);
218 Vmath::Vvtvp(nq, inarray[vz], 1, m_rotMat[2], 1, outarray[vx],
219 1, outarray[vx], 1);
220 Vmath::Vvtvvtp(nq, inarray[vx], 1, m_rotMat[3], 1, inarray[vy],
221 1, m_rotMat[4], 1, outarray[vy], 1);
222 Vmath::Vvtvp(nq, inarray[vz], 1, m_rotMat[5], 1, outarray[vy],
223 1, outarray[vy], 1);
224 Vmath::Vvtvvtp(nq, inarray[vx], 1, m_rotMat[6], 1, inarray[vy],
225 1, m_rotMat[7], 1, outarray[vz], 1);
226 Vmath::Vvtvp(nq, inarray[vz], 1, m_rotMat[8], 1, outarray[vz],
227 1, outarray[vz], 1);
228 break;
229 }
230
231 default:
232 ASSERTL1(false, "Invalid space dimension.");
233 break;
234 }
235 }
236}
237
238/**
239 * @brief Rotate a vector field from trace normal.
240 *
241 * This function performs a rotation of the triad of vector components
242 * provided in inarray so that the first component aligns with the
243 * Cartesian components; it performs the inverse operation of
244 * RiemannSolver::rotateToNormal.
245 */
247 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
248 const Array<OneD, const Array<OneD, NekDouble>> &normals,
249 const Array<OneD, const Array<OneD, NekDouble>> &vecLocs,
251{
252 for (int i = 0; i < inarray.size(); ++i)
253 {
254 Vmath::Vcopy(inarray[i].size(), inarray[i], 1, outarray[i], 1);
255 }
256
257 for (int i = 0; i < vecLocs.size(); i++)
258 {
259 ASSERTL1(vecLocs[i].size() == normals.size(),
260 "vecLocs[i] element count mismatch");
261
262 switch (normals.size())
263 {
264 case 1:
265 { // do nothing
266 const int nq = normals[0].size();
267 const int vx = (int)vecLocs[i][0];
268 Vmath::Vmul(nq, inarray[vx], 1, normals[0], 1, outarray[vx], 1);
269 break;
270 }
271 case 2:
272 {
273 const int nq = normals[0].size();
274 const int vx = (int)vecLocs[i][0];
275 const int vy = (int)vecLocs[i][1];
276
277 Vmath::Vmul(nq, inarray[vy], 1, normals[1], 1, outarray[vx], 1);
278 Vmath::Vvtvm(nq, inarray[vx], 1, normals[0], 1, outarray[vx], 1,
279 outarray[vx], 1);
280 Vmath::Vmul(nq, inarray[vx], 1, normals[1], 1, outarray[vy], 1);
281 Vmath::Vvtvp(nq, inarray[vy], 1, normals[0], 1, outarray[vy], 1,
282 outarray[vy], 1);
283 break;
284 }
285
286 case 3:
287 {
288 const int nq = normals[0].size();
289 const int vx = (int)vecLocs[i][0];
290 const int vy = (int)vecLocs[i][1];
291 const int vz = (int)vecLocs[i][2];
292
293 Vmath::Vvtvvtp(nq, inarray[vx], 1, m_rotMat[0], 1, inarray[vy],
294 1, m_rotMat[3], 1, outarray[vx], 1);
295 Vmath::Vvtvp(nq, inarray[vz], 1, m_rotMat[6], 1, outarray[vx],
296 1, outarray[vx], 1);
297 Vmath::Vvtvvtp(nq, inarray[vx], 1, m_rotMat[1], 1, inarray[vy],
298 1, m_rotMat[4], 1, outarray[vy], 1);
299 Vmath::Vvtvp(nq, inarray[vz], 1, m_rotMat[7], 1, outarray[vy],
300 1, outarray[vy], 1);
301 Vmath::Vvtvvtp(nq, inarray[vx], 1, m_rotMat[2], 1, inarray[vy],
302 1, m_rotMat[5], 1, outarray[vz], 1);
303 Vmath::Vvtvp(nq, inarray[vz], 1, m_rotMat[8], 1, outarray[vz],
304 1, outarray[vz], 1);
305 break;
306 }
307
308 default:
309 ASSERTL1(false, "Invalid space dimension.");
310 break;
311 }
312 }
313}
314
315/**
316 * @brief Determine whether a scalar has been defined in #m_scalars.
317 *
318 * @param name Scalar name.
319 */
321{
322 return m_scalars.find(name) != m_scalars.end();
323}
324
325/**
326 * @brief Determine whether a vector has been defined in #m_vectors.
327 *
328 * @param name Vector name.
329 */
331{
332 return m_vectors.find(name) != m_vectors.end();
333}
334
335/**
336 * @brief Determine whether a parameter has been defined in #m_params.
337 *
338 * @param name Parameter name.
339 */
341{
342 return m_params.find(name) != m_params.end();
343}
344
345/**
346 * @brief Determine whether a scalar has been defined in #m_auxScal.
347 *
348 * @param name Scalar name.
349 */
351{
352 return m_auxScal.find(name) != m_auxScal.end();
353}
354
355/**
356 * @brief Determine whether a vector has been defined in #m_auxVec.
357 *
358 * @param name Vector name.
359 */
361{
362 return m_auxVec.find(name) != m_auxVec.end();
363}
364
365/**
366 * @brief Generate rotation matrices for 3D expansions.
367 */
369 const Array<OneD, const Array<OneD, NekDouble>> &normals)
370{
371 Array<OneD, NekDouble> xdir(3, 0.0);
373 NekDouble tmp[9];
374 const int nq = normals[0].size();
375 int i, j;
376 xdir[0] = 1.0;
377
378 // Allocate storage for rotation matrices.
380
381 for (i = 0; i < 9; ++i)
382 {
384 }
385 for (i = 0; i < normals[0].size(); ++i)
386 {
387 // Generate matrix which takes us from (1,0,0) vector to trace
388 // normal.
389 tn[0] = normals[0][i];
390 tn[1] = normals[1][i];
391 tn[2] = normals[2][i];
392 FromToRotation(tn, xdir, tmp);
393
394 for (j = 0; j < 9; ++j)
395 {
396 m_rotMat[j][i] = tmp[j];
397 }
398 }
399}
400
401/**
402 * @brief A function for creating a rotation matrix that rotates a
403 * vector @a from into another vector @a to.
404 *
405 * Authors: Tomas Möller, John Hughes
406 * "Efficiently Building a Matrix to Rotate One Vector to
407 * Another" Journal of Graphics Tools, 4(4):1-4, 1999
408 *
409 * @param from Normalised 3-vector to rotate from.
410 * @param to Normalised 3-vector to rotate to.
411 * @param out Resulting 3x3 rotation matrix (row-major order).
412 */
415 NekDouble *mat)
416{
417 NekDouble v[3];
418 NekDouble e, h, f;
419
420 CROSS(v, from, to);
421 e = DOT(from, to);
422 f = (e < 0) ? -e : e;
423 if (f > 1.0 - EPSILON)
424 {
425 NekDouble u[3], v[3];
426 NekDouble x[3];
427 NekDouble c1, c2, c3;
428 int i, j;
429
430 x[0] = (from[0] > 0.0) ? from[0] : -from[0];
431 x[1] = (from[1] > 0.0) ? from[1] : -from[1];
432 x[2] = (from[2] > 0.0) ? from[2] : -from[2];
433
434 if (x[0] < x[1])
435 {
436 if (x[0] < x[2])
437 {
438 x[0] = 1.0;
439 x[1] = x[2] = 0.0;
440 }
441 else
442 {
443 x[2] = 1.0;
444 x[0] = x[1] = 0.0;
445 }
446 }
447 else
448 {
449 if (x[1] < x[2])
450 {
451 x[1] = 1.0;
452 x[0] = x[2] = 0.0;
453 }
454 else
455 {
456 x[2] = 1.0;
457 x[0] = x[1] = 0.0;
458 }
459 }
460
461 u[0] = x[0] - from[0];
462 u[1] = x[1] - from[1];
463 u[2] = x[2] - from[2];
464 v[0] = x[0] - to[0];
465 v[1] = x[1] - to[1];
466 v[2] = x[2] - to[2];
467
468 c1 = 2.0 / DOT(u, u);
469 c2 = 2.0 / DOT(v, v);
470 c3 = c1 * c2 * DOT(u, v);
471
472 for (i = 0; i < 3; i++)
473 {
474 for (j = 0; j < 3; j++)
475 {
476 mat[3 * i + j] =
477 -c1 * u[i] * u[j] - c2 * v[i] * v[j] + c3 * v[i] * u[j];
478 }
479 mat[i + 3 * i] += 1.0;
480 }
481 }
482 else
483 {
484 NekDouble hvx, hvz, hvxy, hvxz, hvyz;
485 h = 1.0 / (1.0 + e);
486 hvx = h * v[0];
487 hvz = h * v[2];
488 hvxy = hvx * v[1];
489 hvxz = hvx * v[2];
490 hvyz = hvz * v[1];
491 mat[0] = e + hvx * v[0];
492 mat[1] = hvxy - v[2];
493 mat[2] = hvxz + v[1];
494 mat[3] = hvxy + v[2];
495 mat[4] = e + h * v[1] * v[1];
496 mat[5] = hvyz - v[0];
497 mat[6] = hvxz - v[1];
498 mat[7] = hvyz + v[0];
499 mat[8] = e + hvz * v[2];
500 }
501}
502
503/**
504 * @brief Calculate the flux jacobian of Fwd and Bwd
505 *
506 * @param Fwd Forwards trace space.
507 * @param Bwd Backwards trace space.
508 * @param flux Resultant flux along trace space.
509 */
511 const int nDim, const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
512 const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
514{
515 int nPts = Fwd[0].size();
516
518 {
519 ASSERTL1(CheckVectors("N"), "N not defined.");
520 ASSERTL1(CheckAuxVec("vecLocs"), "vecLocs not defined.");
522 m_vectors["N"]();
524 m_auxVec["vecLocs"]();
525
526 v_CalcFluxJacobian(nDim, Fwd, Bwd, normals, FJac, BJac);
527 }
528 else
529 {
531 for (int i = 0; i < nDim; i++)
532 {
533 normals[i] = Array<OneD, NekDouble>(nPts, 0.0);
534 }
535 Vmath::Fill(nPts, 1.0, normals[0], 1);
536
537 v_CalcFluxJacobian(nDim, Fwd, Bwd, normals, FJac, BJac);
538 }
539}
540
542 [[maybe_unused]] const int nDim,
543 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
544 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
545 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>> &normals,
546 [[maybe_unused]] DNekBlkMatSharedPtr &FJac,
547 [[maybe_unused]] DNekBlkMatSharedPtr &BJac)
548{
549 NEKERROR(ErrorUtil::efatal, "v_CalcFluxJacobian not specified.");
550}
551
552} // 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.
Definition: NekFactory.hpp:104
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()
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