Nektar++
StdExpansion2D.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: StdExpansion2D.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: Daughter of StdExpansion. This class contains routine
32 // which are common to 2D expansion. Typically this inolves physiocal
33 // space operations.
34 //
35 ///////////////////////////////////////////////////////////////////////////////
36 
38 
39 #ifdef max
40 #undef max
41 #endif
42 
43 namespace Nektar
44 {
45 namespace StdRegions
46 {
48 {
49 }
50 
52  const LibUtilities::BasisKey &Bb)
53  : StdExpansion(numcoeffs, 2, Ba, Bb)
54 {
55 }
56 
58 {
59 }
60 
62 {
63 }
64 
65 //----------------------------
66 // Differentiation Methods
67 //----------------------------
69  const Array<OneD, const NekDouble> &inarray,
70  Array<OneD, NekDouble> &outarray_d0, Array<OneD, NekDouble> &outarray_d1)
71 {
72  int nquad0 = m_base[0]->GetNumPoints();
73  int nquad1 = m_base[1]->GetNumPoints();
74 
75  if (outarray_d0.size() > 0) // calculate du/dx_0
76  {
77  DNekMatSharedPtr D0 = m_base[0]->GetD();
78  if (inarray.data() == outarray_d0.data())
79  {
80  Array<OneD, NekDouble> wsp(nquad0 * nquad1);
81  Vmath::Vcopy(nquad0 * nquad1, inarray.get(), 1, wsp.get(), 1);
82  Blas::Dgemm('N', 'N', nquad0, nquad1, nquad0, 1.0,
83  &(D0->GetPtr())[0], nquad0, &wsp[0], nquad0, 0.0,
84  &outarray_d0[0], nquad0);
85  }
86  else
87  {
88  Blas::Dgemm('N', 'N', nquad0, nquad1, nquad0, 1.0,
89  &(D0->GetPtr())[0], nquad0, &inarray[0], nquad0, 0.0,
90  &outarray_d0[0], nquad0);
91  }
92  }
93 
94  if (outarray_d1.size() > 0) // calculate du/dx_1
95  {
96  DNekMatSharedPtr D1 = m_base[1]->GetD();
97  if (inarray.data() == outarray_d1.data())
98  {
99  Array<OneD, NekDouble> wsp(nquad0 * nquad1);
100  Vmath::Vcopy(nquad0 * nquad1, inarray.get(), 1, wsp.get(), 1);
101  Blas::Dgemm('N', 'T', nquad0, nquad1, nquad1, 1.0, &wsp[0], nquad0,
102  &(D1->GetPtr())[0], nquad1, 0.0, &outarray_d1[0],
103  nquad0);
104  }
105  else
106  {
107  Blas::Dgemm('N', 'T', nquad0, nquad1, nquad1, 1.0, &inarray[0],
108  nquad0, &(D1->GetPtr())[0], nquad1, 0.0,
109  &outarray_d1[0], nquad0);
110  }
111  }
112 }
113 
115  const Array<OneD, const NekDouble> &coords,
116  const Array<OneD, const NekDouble> &physvals)
117 {
118  ASSERTL2(coords[0] > -1 - NekConstants::kNekZeroTol, "coord[0] < -1");
119  ASSERTL2(coords[0] < 1 + NekConstants::kNekZeroTol, "coord[0] > 1");
120  ASSERTL2(coords[1] > -1 - NekConstants::kNekZeroTol, "coord[1] < -1");
121  ASSERTL2(coords[1] < 1 + NekConstants::kNekZeroTol, "coord[1] > 1");
122 
123  Array<OneD, NekDouble> coll(2);
124  LocCoordToLocCollapsed(coords, coll);
125 
126  const int nq0 = m_base[0]->GetNumPoints();
127  const int nq1 = m_base[1]->GetNumPoints();
128 
129  Array<OneD, NekDouble> wsp(nq1);
130  for (int i = 0; i < nq1; ++i)
131  {
132  wsp[i] = StdExpansion::BaryEvaluate<0>(coll[0], &physvals[0] + i * nq0);
133  }
134 
135  return StdExpansion::BaryEvaluate<1>(coll[1], &wsp[0]);
136 }
137 
140  const Array<OneD, const NekDouble> &physvals)
141 {
142  NekDouble val;
143  int i;
144  int nq0 = m_base[0]->GetNumPoints();
145  int nq1 = m_base[1]->GetNumPoints();
146  Array<OneD, NekDouble> wsp1(nq1);
147 
148  // interpolate first coordinate direction
149  for (i = 0; i < nq1; ++i)
150  {
151  wsp1[i] =
152  Blas::Ddot(nq0, &(I[0]->GetPtr())[0], 1, &physvals[i * nq0], 1);
153  }
154 
155  // interpolate in second coordinate direction
156  val = Blas::Ddot(nq1, I[1]->GetPtr(), 1, wsp1, 1);
157 
158  return val;
159 }
160 
162  const Array<OneD, NekDouble> &coord,
163  const Array<OneD, const NekDouble> &inarray,
164  std::array<NekDouble, 3> &firstOrderDerivs)
165 {
166  boost::ignore_unused(coord, inarray, firstOrderDerivs);
167  return 0;
168 }
169 
170 //////////////////////////////
171 // Integration Methods
172 //////////////////////////////
173 
177 {
178  int i;
179  NekDouble Int = 0.0;
180  int nquad0 = m_base[0]->GetNumPoints();
181  int nquad1 = m_base[1]->GetNumPoints();
182  Array<OneD, NekDouble> tmp(nquad0 * nquad1);
183 
184  // multiply by integration constants
185  for (i = 0; i < nquad1; ++i)
186  {
187  Vmath::Vmul(nquad0, &inarray[0] + i * nquad0, 1, w0.get(), 1,
188  &tmp[0] + i * nquad0, 1);
189  }
190 
191  for (i = 0; i < nquad0; ++i)
192  {
193  Vmath::Vmul(nquad1, &tmp[0] + i, nquad0, w1.get(), 1, &tmp[0] + i,
194  nquad0);
195  }
196  Int = Vmath::Vsum(nquad0 * nquad1, tmp, 1);
197 
198  return Int;
199 }
200 
202  const Array<OneD, const NekDouble> &base0,
203  const Array<OneD, const NekDouble> &base1,
204  const Array<OneD, const NekDouble> &inarray,
206  bool doCheckCollDir0, bool doCheckCollDir1)
207 {
208  v_BwdTrans_SumFacKernel(base0, base1, inarray, outarray, wsp,
209  doCheckCollDir0, doCheckCollDir1);
210 }
211 
213  const Array<OneD, const NekDouble> &base0,
214  const Array<OneD, const NekDouble> &base1,
215  const Array<OneD, const NekDouble> &inarray,
217  bool doCheckCollDir0, bool doCheckCollDir1)
218 {
219  v_IProductWRTBase_SumFacKernel(base0, base1, inarray, outarray, wsp,
220  doCheckCollDir0, doCheckCollDir1);
221 }
222 
224 {
225  ASSERTL1((dir == 0) || (dir == 1), "Invalid direction.");
226 
227  int nquad0 = m_base[0]->GetNumPoints();
228  int nquad1 = m_base[1]->GetNumPoints();
229  int nqtot = nquad0 * nquad1;
230  int nmodes0 = m_base[0]->GetNumModes();
231 
232  Array<OneD, NekDouble> tmp1(2 * nqtot + m_ncoeffs + nmodes0 * nquad1, 0.0);
233  Array<OneD, NekDouble> tmp3(tmp1 + 2 * nqtot);
234  Array<OneD, NekDouble> tmp4(tmp1 + 2 * nqtot + m_ncoeffs);
235 
236  switch (dir)
237  {
238  case 0:
239  for (int i = 0; i < nqtot; i++)
240  {
241  tmp1[i] = 1.0;
242  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
243  m_base[1]->GetBdata(), tmp1, tmp3,
244  tmp4, false, true);
245  tmp1[i] = 0.0;
246 
247  for (int j = 0; j < m_ncoeffs; j++)
248  {
249  (*mat)(j, i) = tmp3[j];
250  }
251  }
252  break;
253  case 1:
254  for (int i = 0; i < nqtot; i++)
255  {
256  tmp1[i] = 1.0;
257  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
258  m_base[1]->GetDbdata(), tmp1, tmp3,
259  tmp4, true, false);
260  tmp1[i] = 0.0;
261 
262  for (int j = 0; j < m_ncoeffs; j++)
263  {
264  (*mat)(j, i) = tmp3[j];
265  }
266  }
267  break;
268  default:
269  NEKERROR(ErrorUtil::efatal, "Not a 2D expansion.");
270  break;
271  }
272 }
273 
275  const Array<OneD, const NekDouble> &inarray,
276  Array<OneD, NekDouble> &outarray, const StdRegions::StdMatrixKey &mkey)
277 {
278  if (mkey.GetNVarCoeff() == 0 &&
281  {
282  using std::max;
283 
284  // This implementation is only valid when there are no
285  // coefficients associated to the Laplacian operator
286  int nquad0 = m_base[0]->GetNumPoints();
287  int nquad1 = m_base[1]->GetNumPoints();
288  int nqtot = nquad0 * nquad1;
289  int nmodes0 = m_base[0]->GetNumModes();
290  int nmodes1 = m_base[1]->GetNumModes();
291  int wspsize =
292  max(max(max(nqtot, m_ncoeffs), nquad1 * nmodes0), nquad0 * nmodes1);
293 
294  const Array<OneD, const NekDouble> &base0 = m_base[0]->GetBdata();
295  const Array<OneD, const NekDouble> &base1 = m_base[1]->GetBdata();
296 
297  // Allocate temporary storage
298  Array<OneD, NekDouble> wsp0(4 * wspsize); // size wspsize
299  Array<OneD, NekDouble> wsp1(wsp0 + wspsize); // size 3*wspsize
300 
301  if (!(m_base[0]->Collocation() && m_base[1]->Collocation()))
302  {
303  // LAPLACIAN MATRIX OPERATION
304  // wsp0 = u = B * u_hat
305  // wsp1 = du_dxi1 = D_xi1 * wsp0 = D_xi1 * u
306  // wsp2 = du_dxi2 = D_xi2 * wsp0 = D_xi2 * u
307  BwdTrans_SumFacKernel(base0, base1, inarray, wsp0, wsp1, true,
308  true);
309  LaplacianMatrixOp_MatFree_Kernel(wsp0, outarray, wsp1);
310  }
311  else
312  {
313  LaplacianMatrixOp_MatFree_Kernel(inarray, outarray, wsp1);
314  }
315  }
316  else
317  {
319  mkey);
320  }
321 }
322 
324  const Array<OneD, const NekDouble> &inarray,
325  Array<OneD, NekDouble> &outarray, const StdRegions::StdMatrixKey &mkey)
326 {
327  if (mkey.GetNVarCoeff() == 0 &&
330  {
331  using std::max;
332 
333  int nquad0 = m_base[0]->GetNumPoints();
334  int nquad1 = m_base[1]->GetNumPoints();
335  int nqtot = nquad0 * nquad1;
336  int nmodes0 = m_base[0]->GetNumModes();
337  int nmodes1 = m_base[1]->GetNumModes();
338  int wspsize =
339  max(max(max(nqtot, m_ncoeffs), nquad1 * nmodes0), nquad0 * nmodes1);
341 
342  const Array<OneD, const NekDouble> &base0 = m_base[0]->GetBdata();
343  const Array<OneD, const NekDouble> &base1 = m_base[1]->GetBdata();
344 
345  // Allocate temporary storage
346  Array<OneD, NekDouble> wsp0(5 * wspsize); // size wspsize
347  Array<OneD, NekDouble> wsp1(wsp0 + wspsize); // size wspsize
348  Array<OneD, NekDouble> wsp2(wsp0 + 2 * wspsize); // size 3*wspsize
349 
350  if (!(m_base[0]->Collocation() && m_base[1]->Collocation()))
351  {
352  // MASS MATRIX OPERATION
353  // The following is being calculated:
354  // wsp0 = B * u_hat = u
355  // wsp1 = W * wsp0
356  // outarray = B^T * wsp1 = B^T * W * B * u_hat = M * u_hat
357  BwdTrans_SumFacKernel(base0, base1, inarray, wsp0, wsp2, true,
358  true);
359  MultiplyByQuadratureMetric(wsp0, wsp1);
360  IProductWRTBase_SumFacKernel(base0, base1, wsp1, outarray, wsp2,
361  true, true);
362 
363  LaplacianMatrixOp_MatFree_Kernel(wsp0, wsp1, wsp2);
364  }
365  else
366  {
367  MultiplyByQuadratureMetric(inarray, outarray);
368  LaplacianMatrixOp_MatFree_Kernel(inarray, wsp1, wsp2);
369  }
370 
371  // outarray = lambda * outarray + wsp1
372  // = (lambda * M + L ) * u_hat
373  Vmath::Svtvp(m_ncoeffs, lambda, &outarray[0], 1, &wsp1[0], 1,
374  &outarray[0], 1);
375  }
376  else
377  {
379  mkey);
380  }
381 }
382 
383 void StdExpansion2D::v_GetTraceCoeffMap(const unsigned int traceid,
384  Array<OneD, unsigned int> &maparray)
385 {
386  boost::ignore_unused(traceid, maparray);
387 
388  ASSERTL0(false,
389  "This method must be defined at the individual shape level");
390 }
391 
392 /** \brief Determine the mapping to re-orientate the
393  coefficients along the element trace (assumed to align
394  with the standard element) into the orientation of the
395  local trace given by edgeOrient.
396  */
398  const unsigned int eid, Array<OneD, unsigned int> &maparray,
399  Array<OneD, int> &signarray, Orientation edgeOrient, int P, int Q)
400 {
401  // Q is only used in 2D traces.
402  boost::ignore_unused(Q);
403 
404  unsigned int i;
405 
406  int dir;
407  // determine basis direction for edge.
409  {
410  dir = (eid == 0) ? 0 : 1;
411  }
412  else
413  {
414  dir = eid % 2;
415  }
416 
417  int numModes = m_base[dir]->GetNumModes();
418 
419  // P is the desired length of the map
420  P = (P == -1) ? numModes : P;
421 
422  // decalare maparray
423  if (maparray.size() != P)
424  {
425  maparray = Array<OneD, unsigned int>(P);
426  }
427 
428  // fill default mapping as increasing index
429  for (i = 0; i < P; ++i)
430  {
431  maparray[i] = i;
432  }
433 
434  if (signarray.size() != P)
435  {
436  signarray = Array<OneD, int>(P, 1);
437  }
438  else
439  {
440  std::fill(signarray.get(), signarray.get() + P, 1);
441  }
442 
443  // Zero signmap and set maparray to zero if
444  // elemental modes are not as large as trace modes
445  for (i = numModes; i < P; ++i)
446  {
447  signarray[i] = 0.0;
448  maparray[i] = maparray[0];
449  }
450 
451  if (edgeOrient == eBackwards)
452  {
453  const LibUtilities::BasisType bType = GetBasisType(dir);
454 
455  if ((bType == LibUtilities::eModified_A) ||
456  (bType == LibUtilities::eModified_B))
457  {
458  std::swap(maparray[0], maparray[1]);
459 
460  for (i = 3; i < std::min(P, numModes); i += 2)
461  {
462  signarray[i] *= -1;
463  }
464  }
465  else if (bType == LibUtilities::eGLL_Lagrange ||
467  {
468  ASSERTL1(P == numModes, "Different trace space edge dimension "
469  "and element edge dimension not currently "
470  "possible for GLL-Lagrange bases");
471 
472  std::reverse(maparray.get(), maparray.get() + P);
473  }
474  else
475  {
476  ASSERTL0(false, "Mapping not defined for this type of basis");
477  }
478  }
479 }
480 
482  Array<OneD, unsigned int> &maparray,
483  Array<OneD, int> &signarray,
484  Orientation edgeOrient, int P,
485  int Q)
486 {
487  Array<OneD, unsigned int> map1, map2;
488  v_GetTraceCoeffMap(eid, map1);
489  v_GetElmtTraceToTraceMap(eid, map2, signarray, edgeOrient, P, Q);
490 
491  if (maparray.size() != map2.size())
492  {
493  maparray = Array<OneD, unsigned int>(map2.size());
494  }
495 
496  for (int i = 0; i < map2.size(); ++i)
497  {
498  maparray[i] = map1[map2[i]];
499  }
500 }
501 
502 } // namespace StdRegions
503 } // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#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 ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:272
Describes the specification for a Basis.
Definition: Basis.h:50
virtual void v_GenStdMatBwdDeriv(const int dir, DNekMatSharedPtr &mat) override
virtual void v_IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1)=0
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d0, Array< OneD, NekDouble > &outarray_d1)
Calculate the 2D derivative in the local tensor/collapsed coordinate at the physical points.
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
virtual void v_GetTraceToElementMap(const int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient=eForwards, int P=-1, int Q=-1) override
virtual void v_BwdTrans_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1)=0
NekDouble Integral(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &w0, const Array< OneD, const NekDouble > &w1)
void BwdTrans_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
virtual void v_GetTraceCoeffMap(const unsigned int traceid, Array< OneD, unsigned int > &maparray) override
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
void IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
virtual NekDouble v_PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals) override
This function evaluates the expansion at a single (arbitrary) point of the domain.
virtual void v_GetElmtTraceToTraceMap(const unsigned int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient, int P, int Q) override
Determine the mapping to re-orientate the coefficients along the element trace (assumed to align with...
The base class for all shapes.
Definition: StdExpansion.h:71
void LaplacianMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:162
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:373
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:729
void HelmholtzMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void LaplacianMatrixOp_MatFree_Kernel(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
NekDouble GetConstFactor(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:126
bool ConstFactorExists(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:135
static double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = .
Definition: Blas.hpp:182
static void Dgemm(const char &transa, const char &transb, const int &m, const int &n, const int &k, const double &alpha, const double *a, const int &lda, const double *b, const int &ldb, const double &beta, double *c, const int &ldc)
BLAS level 3: Matrix-matrix multiply C = A x B where op(A)[m x k], op(B)[k x n], C[m x n] DGEMM perfo...
Definition: Blas.hpp:368
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:51
@ eGauss_Lagrange
Lagrange Polynomials using the Gauss points.
Definition: BasisType.h:59
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:58
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:50
static const NekDouble kNekZeroTol
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
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 Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
Definition: Vmath.cpp:622
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:895
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255