Nektar++
StdExpansion.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: StdExpansion.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: Definition of methods in class StdExpansion which is
32 // the base class to all expansion shapes
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #include <boost/core/ignore_unused.hpp>
37 
38 #include <LibUtilities/Foundations/ManagerAccess.h> // for BasisManager, etc
40 
41 namespace Nektar
42 {
43 namespace StdRegions
44 {
45 /** \brief Default constructor */
47 {
48 }
49 
50 StdExpansion::StdExpansion(const int numcoeffs, const int numbases,
51  const LibUtilities::BasisKey &Ba,
52  const LibUtilities::BasisKey &Bb,
53  const LibUtilities::BasisKey &Bc)
54  : m_base(numbases), m_elmt_id(0), m_ncoeffs(numcoeffs),
55  m_stdMatrixManager(std::bind(&StdExpansion::CreateStdMatrix, this,
56  std::placeholders::_1),
57  std::string("StdExpansionStdMatrix")),
58  m_stdStaticCondMatrixManager(
59  std::bind(&StdExpansion::CreateStdStaticCondMatrix, this,
60  std::placeholders::_1),
61  std::string("StdExpansionStdStaticCondMatrix"))
62 {
63  switch (m_base.size())
64  {
65  case 3:
67  "NULL Basis attempting to be used.");
69  /* Falls through. */
70  case 2:
72  "NULL Basis attempting to be used.");
74  /* Falls through. */
75  case 1:
77  "NULL Basis attempting to be used.");
79  break;
80  default:
81  break;
82  // ASSERTL0(false, "numbases incorrectly specified");
83  };
84 
85 } // end constructor
86 
88  : std::enable_shared_from_this<StdExpansion>(T), m_base(T.m_base),
89  m_elmt_id(T.m_elmt_id), m_ncoeffs(T.m_ncoeffs),
90  m_stdMatrixManager(T.m_stdMatrixManager),
91  m_stdStaticCondMatrixManager(T.m_stdStaticCondMatrixManager)
92 {
93 }
94 
95 // Destructor
97 {
98 }
99 
101  const Array<OneD, const NekDouble> &sol)
102 {
103  NekDouble val;
104  int ntot = GetTotPoints();
105  Array<OneD, NekDouble> wsp(ntot);
106 
107  if (sol == NullNekDouble1DArray)
108  {
109  Vmath::Vabs(ntot, phys, 1, wsp, 1);
110  }
111  else
112  {
113  Vmath::Vsub(ntot, sol, 1, phys, 1, wsp, 1);
114  Vmath::Vabs(ntot, wsp, 1, wsp, 1);
115  }
116 
117  val = Vmath::Vamax(ntot, wsp, 1);
118  return val;
119 }
120 
122  const Array<OneD, const NekDouble> &sol)
123 {
124  NekDouble val;
125  int ntot = GetTotPoints();
126  Array<OneD, NekDouble> wsp(ntot);
127 
128  if (sol.size() == 0)
129  {
130  Vmath::Vmul(ntot, phys, 1, phys, 1, wsp, 1);
131  }
132  else
133  {
134  Vmath::Vsub(ntot, sol, 1, phys, 1, wsp, 1);
135  Vmath::Vmul(ntot, wsp, 1, wsp, 1, wsp, 1);
136  }
137 
138  val = v_Integral(wsp);
139 
140  // if val too small, sqrt returns nan.
142  {
143  return 0.0;
144  }
145  else
146  {
147  return sqrt(val);
148  }
149 }
150 
152  const Array<OneD, const NekDouble> &sol)
153 {
154  int i;
155  NekDouble val;
156  int ntot = GetTotPoints();
157  int coordim = v_GetCoordim();
158  Array<OneD, NekDouble> wsp(3 * ntot);
159  Array<OneD, NekDouble> wsp_deriv = wsp + ntot;
160  Array<OneD, NekDouble> sum = wsp_deriv + ntot;
161 
162  if (sol == NullNekDouble1DArray)
163  {
164  Vmath::Vcopy(ntot, phys, 1, wsp, 1);
165  Vmath::Vmul(ntot, phys, 1, phys, 1, sum, 1);
166  }
167  else
168  {
169  Vmath::Vsub(ntot, sol, 1, phys, 1, wsp, 1);
170  Vmath::Vmul(ntot, wsp, 1, wsp, 1, sum, 1);
171  }
172 
173  for (i = 0; i < coordim; ++i)
174  {
175  v_PhysDeriv(i, wsp, wsp_deriv);
176  Vmath::Vvtvp(ntot, wsp_deriv, 1, wsp_deriv, 1, sum, 1, sum, 1);
177  }
178 
179  val = sqrt(v_Integral(sum));
180 
181  return val;
182 }
183 
185  const StdMatrixKey &mkey)
186 {
187  DNekBlkMatSharedPtr returnval;
188 
189  DNekMatSharedPtr mat = GetStdMatrix(mkey);
190  int nbdry = NumBndryCoeffs(); // also checks to see if this is a boundary
191  // interior decomposed expansion
192  int nint = m_ncoeffs - nbdry;
198 
199  int i, j;
200 
201  Array<OneD, unsigned int> bmap(nbdry);
202  Array<OneD, unsigned int> imap(nint);
203  GetBoundaryMap(bmap);
204  GetInteriorMap(imap);
205 
206  for (i = 0; i < nbdry; ++i)
207  {
208  for (j = 0; j < nbdry; ++j)
209  {
210  (*A)(i, j) = (*mat)(bmap[i], bmap[j]);
211  }
212 
213  for (j = 0; j < nint; ++j)
214  {
215  (*B)(i, j) = (*mat)(bmap[i], imap[j]);
216  }
217  }
218 
219  for (i = 0; i < nint; ++i)
220  {
221  for (j = 0; j < nbdry; ++j)
222  {
223  (*C)(i, j) = (*mat)(imap[i], bmap[j]);
224  }
225 
226  for (j = 0; j < nint; ++j)
227  {
228  (*D)(i, j) = (*mat)(imap[i], imap[j]);
229  }
230  }
231 
232  // Calculate static condensed system
233  if (nint)
234  {
235  D->Invert();
236  (*B) = (*B) * (*D);
237  (*A) = (*A) - (*B) * (*C);
238  }
239 
240  // set up block matrix system
241  Array<OneD, unsigned int> exp_size(2);
242  exp_size[0] = nbdry;
243  exp_size[1] = nint;
244  returnval =
246 
247  returnval->SetBlock(0, 0, A);
248  returnval->SetBlock(0, 1, B);
249  returnval->SetBlock(1, 0, C);
250  returnval->SetBlock(1, 1, D);
251 
252  return returnval;
253 }
254 
256 {
257  int i;
258  DNekMatSharedPtr returnval;
259 
260  switch (mkey.GetMatrixType())
261  {
262  case eInvMass:
263  {
264  StdMatrixKey masskey(eMass, mkey.GetShapeType(), *this,
266  mkey.GetNodalPointsType());
267  DNekMatSharedPtr mmat = GetStdMatrix(masskey);
268 
270  *mmat); // Populate standard mass matrix.
271  returnval->Invert();
272  }
273  break;
274  case eInvNBasisTrans:
275  {
276  StdMatrixKey tmpkey(eNBasisTrans, mkey.GetShapeType(), *this,
278  mkey.GetNodalPointsType());
279  DNekMatSharedPtr tmpmat = GetStdMatrix(tmpkey);
281  *tmpmat); // Populate matrix.
282  returnval->Invert();
283  }
284  break;
285  case eBwdMat:
286  {
287  int nq = GetTotPoints();
289  Array<OneD, NekDouble> tmpout(nq);
290 
291  returnval =
293  Array<OneD, NekDouble> Bwd_data = returnval->GetPtr();
294 
296  this->DetShapeType(), *this);
297  DNekMatSharedPtr MatBwdTrans = GetStdMatrix(matkey);
298  Array<OneD, NekDouble> BwdTrans_data = MatBwdTrans->GetPtr();
299 
300  for (int i = 0; i < m_ncoeffs; ++i)
301  {
302  Array<OneD, NekDouble> tmpinn = BwdTrans_data + nq * i;
303  Array<OneD, NekDouble> tmpout = Bwd_data + i;
304 
305  Vmath::Vcopy(nq, tmpinn, 1, tmpout, m_ncoeffs);
306  }
307  }
308  break;
309  case eBwdTrans:
310  {
311  int nq = GetTotPoints();
313  Array<OneD, NekDouble> tmpout(nq);
314 
315  returnval =
317 
318  for (int i = 0; i < m_ncoeffs; ++i)
319  {
320  Vmath::Zero(m_ncoeffs, tmpin, 1);
321  tmpin[i] = 1.0;
322 
323  BwdTrans_SumFac(tmpin, tmpout);
324 
325  Vmath::Vcopy(nq, tmpout.get(), 1,
326  returnval->GetRawPtr() + i * nq, 1);
327  }
328  }
329  break;
330  case eIProductWRTBase:
331  {
332  int nq = GetTotPoints();
333  Array<OneD, NekDouble> tmpin(nq);
335 
336  returnval =
338 
339  for (i = 0; i < nq; ++i)
340  {
341  Vmath::Zero(nq, tmpin, 1);
342  tmpin[i] = 1.0;
343 
344  IProductWRTBase_SumFac(tmpin, tmpout);
345 
346  Vmath::Vcopy(m_ncoeffs, tmpout.get(), 1,
347  returnval->GetRawPtr() + i * m_ncoeffs, 1);
348  }
349  }
350  break;
352  {
353  int nq = GetTotPoints();
354  Array<OneD, NekDouble> tmpin(nq);
356 
357  returnval =
359 
360  for (i = 0; i < nq; ++i)
361  {
362  Vmath::Zero(nq, tmpin, 1);
363  tmpin[i] = 1.0;
364 
365  IProductWRTDerivBase_SumFac(0, tmpin, tmpout);
366 
367  Vmath::Vcopy(m_ncoeffs, tmpout.get(), 1,
368  returnval->GetRawPtr() + i * m_ncoeffs, 1);
369  }
370  }
371  break;
373  {
374  int nq = GetTotPoints();
375  Array<OneD, NekDouble> tmpin(nq);
377 
378  returnval =
380 
381  for (i = 0; i < nq; ++i)
382  {
383  Vmath::Zero(nq, tmpin, 1);
384  tmpin[i] = 1.0;
385 
386  IProductWRTDerivBase_SumFac(1, tmpin, tmpout);
387 
388  Vmath::Vcopy(m_ncoeffs, tmpout.get(), 1,
389  returnval->GetRawPtr() + i * m_ncoeffs, 1);
390  }
391  }
392  break;
394  {
395  int nq = GetTotPoints();
396  Array<OneD, NekDouble> tmpin(nq);
398 
399  returnval =
401 
402  for (i = 0; i < nq; ++i)
403  {
404  Vmath::Zero(nq, tmpin, 1);
405  tmpin[i] = 1.0;
406 
407  IProductWRTDerivBase_SumFac(2, tmpin, tmpout);
408 
409  Vmath::Vcopy(m_ncoeffs, tmpout.get(), 1,
410  returnval->GetRawPtr() + i * m_ncoeffs, 1);
411  }
412  }
413  break;
414  case eDerivBase0:
415  {
416  int nq = GetTotPoints();
417  returnval =
419  GenStdMatBwdDeriv(0, returnval);
420  }
421  break;
422  case eDerivBase1:
423  {
424  int nq = GetTotPoints();
425  returnval =
427  GenStdMatBwdDeriv(1, returnval);
428  }
429  break;
430  case eDerivBase2:
431  {
432  int nq = GetTotPoints();
433  returnval =
435  GenStdMatBwdDeriv(2, returnval);
436  }
437  break;
438  case eEquiSpacedToCoeffs:
439  {
440  // check to see if equispaced basis
441  int nummodes = m_base[0]->GetNumModes();
442  bool equispaced = true;
443  for (int i = 1; i < m_base.size(); ++i)
444  {
445  if (m_base[i]->GetNumModes() != nummodes)
446  {
447  equispaced = false;
448  }
449  }
450 
451  ASSERTL0(equispaced,
452  "Currently need to have same num modes in all "
453  "directionmodes to use EquiSpacedToCoeff method");
454 
455  int ntot = GetTotPoints();
456  Array<OneD, NekDouble> qmode(ntot);
458 
459  returnval =
461  for (int i = 0; i < m_ncoeffs; ++i)
462  {
463  // Get mode at quadrature points
464  FillMode(i, qmode);
465 
466  // interpolate to equi spaced
467  PhysInterpToSimplexEquiSpaced(qmode, emode, nummodes);
468 
469  // fill matrix
470  Vmath::Vcopy(m_ncoeffs, &emode[0], 1,
471  returnval->GetRawPtr() + i * m_ncoeffs, 1);
472  }
473  // invert matrix
474  returnval->Invert();
475  }
476  break;
477  case eMass:
478  case eHelmholtz:
479  case eLaplacian:
480  case eLaplacian00:
481  case eLaplacian01:
482  case eLaplacian02:
483  case eLaplacian11:
484  case eLaplacian12:
485  case eLaplacian22:
486  case eWeakDeriv0:
487  case eWeakDeriv1:
488  case eWeakDeriv2:
490  case eMassLevelCurvature:
493  {
495  returnval =
497  DNekMat &Mat = *returnval;
498 
499  for (i = 0; i < m_ncoeffs; ++i)
500  {
501  Vmath::Zero(m_ncoeffs, tmp, 1);
502  tmp[i] = 1.0;
503 
504  GeneralMatrixOp_MatFree(tmp, tmp, mkey);
505 
506  Vmath::Vcopy(m_ncoeffs, &tmp[0], 1,
507  &(Mat.GetPtr())[0] + i * m_ncoeffs, 1);
508  }
509  }
510  break;
511  default:
512  {
514  "This type of matrix, " +
515  static_cast<std::string>(
516  MatrixTypeMap[mkey.GetMatrixType()]) +
517  ", can not be created using a general approach");
518  }
519  break;
520  }
521 
522  return returnval;
523 }
524 
526  Array<OneD, NekDouble> &outarray,
527  const StdMatrixKey &mkey)
528 {
529  switch (mkey.GetMatrixType())
530  {
531  case eMass:
532  MassMatrixOp(inarray, outarray, mkey);
533  break;
534  case eWeakDeriv0:
535  WeakDerivMatrixOp(0, inarray, outarray, mkey);
536  break;
537  case eWeakDeriv1:
538  WeakDerivMatrixOp(1, inarray, outarray, mkey);
539  break;
540  case eWeakDeriv2:
541  WeakDerivMatrixOp(2, inarray, outarray, mkey);
542  break;
544  WeakDirectionalDerivMatrixOp(inarray, outarray, mkey);
545  break;
546  case eMassLevelCurvature:
547  MassLevelCurvatureMatrixOp(inarray, outarray, mkey);
548  break;
550  LinearAdvectionDiffusionReactionMatrixOp(inarray, outarray, mkey,
551  false);
552  break;
554  LinearAdvectionDiffusionReactionMatrixOp(inarray, outarray, mkey);
555  break;
556  case eLaplacian:
557  LaplacianMatrixOp(inarray, outarray, mkey);
558  break;
559  case eLaplacian00:
560  LaplacianMatrixOp(0, 0, inarray, outarray, mkey);
561  break;
562  case eLaplacian01:
563  LaplacianMatrixOp(0, 1, inarray, outarray, mkey);
564  break;
565  case eLaplacian02:
566  LaplacianMatrixOp(0, 2, inarray, outarray, mkey);
567  break;
568  case eLaplacian10:
569  LaplacianMatrixOp(1, 0, inarray, outarray, mkey);
570  break;
571  case eLaplacian11:
572  LaplacianMatrixOp(1, 1, inarray, outarray, mkey);
573  break;
574  case eLaplacian12:
575  LaplacianMatrixOp(1, 2, inarray, outarray, mkey);
576  break;
577  case eLaplacian20:
578  LaplacianMatrixOp(2, 0, inarray, outarray, mkey);
579  break;
580  case eLaplacian21:
581  LaplacianMatrixOp(2, 1, inarray, outarray, mkey);
582  break;
583  case eLaplacian22:
584  LaplacianMatrixOp(2, 2, inarray, outarray, mkey);
585  break;
586  case eHelmholtz:
587  HelmholtzMatrixOp(inarray, outarray, mkey);
588  break;
589  default:
591  "This matrix does not have an operator");
592  break;
593  }
594 }
595 
597  const Array<OneD, const NekDouble> &inarray,
598  Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
599 {
600  switch (mkey.GetMatrixType())
601  {
602  case eMass:
603  MassMatrixOp_MatFree(inarray, outarray, mkey);
604  break;
605  case eWeakDeriv0:
606  WeakDerivMatrixOp_MatFree(0, inarray, outarray, mkey);
607  break;
608  case eWeakDeriv1:
609  WeakDerivMatrixOp_MatFree(1, inarray, outarray, mkey);
610  break;
611  case eWeakDeriv2:
612  WeakDerivMatrixOp_MatFree(2, inarray, outarray, mkey);
613  break;
615  WeakDirectionalDerivMatrixOp_MatFree(inarray, outarray, mkey);
616  break;
617  case eMassLevelCurvature:
618  MassLevelCurvatureMatrixOp_MatFree(inarray, outarray, mkey);
619  break;
622  mkey, false);
623  break;
626  mkey);
627  break;
628  case eLaplacian:
629  LaplacianMatrixOp_MatFree(inarray, outarray, mkey);
630  break;
631  case eLaplacian00:
632  LaplacianMatrixOp_MatFree(0, 0, inarray, outarray, mkey);
633  break;
634  case eLaplacian01:
635  LaplacianMatrixOp_MatFree(0, 1, inarray, outarray, mkey);
636  break;
637  case eLaplacian02:
638  LaplacianMatrixOp_MatFree(0, 2, inarray, outarray, mkey);
639  break;
640  case eLaplacian10:
641  LaplacianMatrixOp_MatFree(1, 0, inarray, outarray, mkey);
642  break;
643  case eLaplacian11:
644  LaplacianMatrixOp_MatFree(1, 1, inarray, outarray, mkey);
645  break;
646  case eLaplacian12:
647  LaplacianMatrixOp_MatFree(1, 2, inarray, outarray, mkey);
648  break;
649  case eLaplacian20:
650  LaplacianMatrixOp_MatFree(2, 0, inarray, outarray, mkey);
651  break;
652  case eLaplacian21:
653  LaplacianMatrixOp_MatFree(2, 1, inarray, outarray, mkey);
654  break;
655  case eLaplacian22:
656  LaplacianMatrixOp_MatFree(2, 2, inarray, outarray, mkey);
657  break;
658  case eHelmholtz:
659  HelmholtzMatrixOp_MatFree(inarray, outarray, mkey);
660  break;
661  default:
663  "This matrix does not have an operator");
664  break;
665  }
666 }
667 
669  const Array<OneD, const NekDouble> &inarray,
670  Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
671 {
672  int nq = GetTotPoints();
673  Array<OneD, NekDouble> tmp(nq);
674 
675  v_BwdTrans(inarray, tmp);
676 
677  if (mkey.HasVarCoeff(eVarCoeffMass))
678  {
679  Vmath::Vmul(nq, mkey.GetVarCoeff(eVarCoeffMass), 1, tmp, 1, tmp, 1);
680  }
681 
682  v_IProductWRTBase(tmp, outarray);
683 }
684 
686  const int k1, const int k2, const Array<OneD, const NekDouble> &inarray,
687  Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
688 {
689  ASSERTL1(k1 >= 0 && k1 < GetCoordim(), "invalid first argument");
690  ASSERTL1(k2 >= 0 && k2 < GetCoordim(), "invalid second argument");
691 
692  int nq = GetTotPoints();
693  Array<OneD, NekDouble> tmp(nq);
694  Array<OneD, NekDouble> dtmp(nq);
695  VarCoeffType varcoefftypes[3][3] = {
699 
700  ConstFactorType constcoefftypes[3][3] = {
704 
705  v_BwdTrans(inarray, tmp);
706  v_PhysDeriv(k2, tmp, dtmp);
707  if (mkey.GetNVarCoeff() && (!mkey.ConstFactorExists(eFactorSVVDiffCoeff)))
708  {
709  if (k1 == k2)
710  {
711  // By default, k1 == k2 has \sigma = 1 (diagonal entries)
712  if (mkey.HasVarCoeff(varcoefftypes[k1][k1]))
713  {
714  Vmath::Vmul(nq, mkey.GetVarCoeff(varcoefftypes[k1][k1]), 1,
715  dtmp, 1, dtmp, 1);
716  }
717  v_IProductWRTDerivBase_SumFac(k1, dtmp, outarray);
718  }
719  else
720  {
721  // By default, k1 != k2 has \sigma = 0 (off-diagonal entries)
722  if (mkey.HasVarCoeff(varcoefftypes[k1][k2]))
723  {
724  Vmath::Vmul(nq, mkey.GetVarCoeff(varcoefftypes[k1][k2]), 1,
725  dtmp, 1, dtmp, 1);
726  v_IProductWRTDerivBase_SumFac(k1, dtmp, outarray);
727  }
728  else
729  {
730  Vmath::Zero(GetNcoeffs(), outarray, 1);
731  }
732  }
733  }
734  else if (mkey.ConstFactorExists(eFactorCoeffD00) &&
736  {
737  if (k1 == k2)
738  {
739  // By default, k1 == k2 has \sigma = 1 (diagonal entries)
740  if (mkey.ConstFactorExists(constcoefftypes[k1][k1]))
741  {
742  Vmath::Smul(nq, mkey.GetConstFactor(constcoefftypes[k1][k1]),
743  dtmp, 1, dtmp, 1);
744  }
745  v_IProductWRTDerivBase(k1, dtmp, outarray);
746  }
747  else
748  {
749  // By default, k1 != k2 has \sigma = 0 (off-diagonal entries)
750  if (mkey.ConstFactorExists(constcoefftypes[k1][k2]))
751  {
752  Vmath::Smul(nq, mkey.GetConstFactor(constcoefftypes[k1][k2]),
753  dtmp, 1, dtmp, 1);
754  v_IProductWRTDerivBase(k1, dtmp, outarray);
755  }
756  else
757  {
758  Vmath::Zero(GetNcoeffs(), outarray, 1);
759  }
760  }
761  }
762  else
763  {
764  // Multiply by svv tensor
766  {
767  Vmath::Vcopy(nq, dtmp, 1, tmp, 1);
768  SVVLaplacianFilter(dtmp, mkey);
769  Vmath::Vadd(nq, tmp, 1, dtmp, 1, dtmp, 1);
770  }
771  v_IProductWRTDerivBase(k1, dtmp, outarray);
772  }
773 }
774 
776  const Array<OneD, const NekDouble> &inarray,
777  Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
778 {
779  const int dim = GetCoordim();
780 
781  int i, j;
782 
784  Array<OneD, NekDouble> store2(m_ncoeffs, 0.0);
785 
786  if ((mkey.GetNVarCoeff() == 0 &&
789  {
790  // just call diagonal matrix form of laplcian operator
791  for (i = 0; i < dim; ++i)
792  {
793  LaplacianMatrixOp(i, i, inarray, store, mkey);
794  Vmath::Vadd(m_ncoeffs, store, 1, store2, 1, store2, 1);
795  }
796  }
797  else
798  {
799  const MatrixType mtype[3][3] = {
803  StdMatrixKeySharedPtr mkeyij;
804 
805  for (i = 0; i < dim; i++)
806  {
807  for (j = 0; j < dim; j++)
808  {
810  mkey, mtype[i][j]);
811  LaplacianMatrixOp(i, j, inarray, store, *mkeyij);
812  Vmath::Vadd(m_ncoeffs, store, 1, store2, 1, store2, 1);
813  }
814  }
815  }
816 
817  Vmath::Vcopy(m_ncoeffs, store2.get(), 1, outarray.get(), 1);
818 }
819 
821  const int k1, const Array<OneD, const NekDouble> &inarray,
822  Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
823 {
825  int nq = GetTotPoints();
826 
827  v_BwdTrans(inarray, tmp);
828  v_PhysDeriv(k1, tmp, tmp);
829 
831  if (mkey.HasVarCoeff(keys[k1]))
832  {
833  Vmath::Vmul(nq, &(mkey.GetVarCoeff(keys[k1]))[0], 1, &tmp[0], 1,
834  &tmp[0], 1);
835  }
836 
837  v_IProductWRTBase(tmp, outarray);
838 }
839 
841  const Array<OneD, const NekDouble> &inarray,
842  Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
843 {
844  int nq = GetTotPoints();
845 
846  Array<OneD, NekDouble> tmp(nq), Dtmp(nq);
847  Array<OneD, NekDouble> Mtmp(nq), Mout(m_ncoeffs);
848 
849  v_BwdTrans(inarray, tmp);
851 
852  v_IProductWRTBase(Dtmp, outarray);
853 
854  // Compte M_{div tv}
855  Vmath::Vmul(nq, &(mkey.GetVarCoeff(eVarCoeffMFDiv))[0], 1, &tmp[0], 1,
856  &Mtmp[0], 1);
857 
858  v_IProductWRTBase(Mtmp, Mout);
859 
860  // Add D_tv + M_{div tv}
861  Vmath::Vadd(m_ncoeffs, &Mout[0], 1, &outarray[0], 1, &outarray[0], 1);
862 }
863 
865  const Array<OneD, const NekDouble> &inarray,
866  Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
867 {
868  boost::ignore_unused(inarray, outarray, mkey);
869 }
870 
872  const Array<OneD, const NekDouble> &inarray,
873  Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey,
874  bool addDiffusionTerm)
875 {
876 
877  int i;
878  int ndir =
879  mkey.GetNVarCoeff(); // assume num.r consts corresponds to directions
880  ASSERTL0(ndir, "Must define at least one advection velocity");
881 
882  NekDouble lambda = mkey.GetConstFactor(eFactorLambda);
883  int totpts = GetTotPoints();
884  Array<OneD, NekDouble> tmp(3 * totpts);
885  Array<OneD, NekDouble> tmp_deriv = tmp + totpts;
886  Array<OneD, NekDouble> tmp_adv = tmp_deriv + totpts;
887 
888  ASSERTL1(ndir <= GetCoordim(),
889  "Number of constants is larger than coordinate dimensions");
890 
891  v_BwdTrans(inarray, tmp);
892 
893  VarCoeffType varcoefftypes[] = {eVarCoeffVelX, eVarCoeffVelY,
894  eVarCoeffVelZ};
895 
896  // calculate u dx + v dy + ..
897  Vmath::Zero(totpts, tmp_adv, 1);
898  for (i = 0; i < ndir; ++i)
899  {
900  v_PhysDeriv(i, tmp, tmp_deriv);
901  Vmath::Vvtvp(totpts, mkey.GetVarCoeff(varcoefftypes[i]), 1, tmp_deriv,
902  1, tmp_adv, 1, tmp_adv, 1);
903  }
904 
905  if (lambda) // add -lambda*u
906  {
907  Vmath::Svtvp(totpts, -lambda, tmp, 1, tmp_adv, 1, tmp_adv, 1);
908  }
909 
910  if (addDiffusionTerm)
911  {
913  StdMatrixKey mkeylap(eLaplacian, DetShapeType(), *this,
914  mkey.GetConstFactors(), mkey.GetVarCoeffs(),
915  mkey.GetNodalPointsType());
916  LaplacianMatrixOp(inarray, lap, mkeylap);
917 
918  v_IProductWRTBase(tmp_adv, outarray);
919  // Lap v - u.grad v + lambda*u
920  // => (grad u, grad v) + u.grad v - lambda*u
921  Vmath::Vadd(m_ncoeffs, lap, 1, outarray, 1, outarray, 1);
922  }
923  else
924  {
925  v_IProductWRTBase(tmp_adv, outarray);
926  }
927 }
928 
930  const Array<OneD, const NekDouble> &inarray,
931  Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
932 {
933  NekDouble lambda = mkey.GetConstFactor(eFactorLambda);
935  StdMatrixKey mkeymass(eMass, DetShapeType(), *this);
936  StdMatrixKey mkeylap(eLaplacian, DetShapeType(), *this,
937  mkey.GetConstFactors(), mkey.GetVarCoeffs(),
938  mkey.GetNodalPointsType());
939 
940  MassMatrixOp(inarray, tmp, mkeymass);
941  LaplacianMatrixOp(inarray, outarray, mkeylap);
942 
943  Blas::Daxpy(m_ncoeffs, lambda, tmp, 1, outarray, 1);
944 }
945 
946 // VIRTUAL INLINE FUNCTIONS FROM HEADER FILE
948  const Array<OneD, const NekDouble> &Lcoord,
949  const Array<OneD, const NekDouble> &physvals)
950 {
951  return v_StdPhysEvaluate(Lcoord, physvals);
952 }
953 
955  const std::vector<unsigned int> &nummodes, int &modes_offset)
956 {
957  boost::ignore_unused(nummodes, modes_offset);
958  NEKERROR(ErrorUtil::efatal, "This function is not defined for this class");
959  return 0;
960 }
961 
964 {
965  boost::ignore_unused(Fx, outarray);
966  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
967 }
968 
972 {
973  boost::ignore_unused(Fx, Fy, outarray);
974  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
975 }
976 
981 {
982  boost::ignore_unused(Fx, Fy, Fz, outarray);
983  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
984 }
985 
987  const Array<OneD, const Array<OneD, NekDouble>> &Fvec,
988  Array<OneD, NekDouble> &outarray)
989 {
990  boost::ignore_unused(Fvec, outarray);
991  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
992 }
993 
995  const LocalRegions::MatrixKey &mkey)
996 {
997  boost::ignore_unused(mkey);
998  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
1000 }
1001 
1003  const LocalRegions::MatrixKey &mkey)
1004 {
1005  boost::ignore_unused(mkey);
1006  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
1007 }
1008 
1011  Array<OneD, NekDouble> &outarray)
1012 {
1013  boost::ignore_unused(dir, inarray, outarray);
1014  NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1015 }
1016 
1019 {
1020  boost::ignore_unused(coeffs, dir);
1021  NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1022 }
1023 
1025  const Array<OneD, const NekDouble> &Lcoord,
1026  const Array<OneD, const NekDouble> &physvals)
1027 
1028 {
1029  boost::ignore_unused(Lcoord, physvals);
1030  NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1031  return 0;
1032 }
1033 
1036 {
1037  boost::ignore_unused(xi, eta);
1038  NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1039 }
1040 
1043 {
1044  boost::ignore_unused(eta, xi);
1045  NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1046 }
1047 
1049  const int k) const
1050 {
1051  boost::ignore_unused(i, k);
1052  ASSERTL0(false, "This function is not valid or not defined");
1054 }
1055 
1057  const int j) const
1058 {
1059  boost::ignore_unused(i, j);
1060  ASSERTL0(false, "This function is not valid or not defined");
1062 }
1063 
1065 {
1066  ASSERTL0(false, "This function is not valid or not defined");
1067 
1069 }
1070 
1071 std::shared_ptr<StdExpansion> StdExpansion::v_GetStdExp(void) const
1072 {
1073  ASSERTL0(false, "This method is not defined for this expansion");
1074  StdExpansionSharedPtr returnval;
1075  return returnval;
1076 }
1077 
1078 std::shared_ptr<StdExpansion> StdExpansion::v_GetLinStdExp(void) const
1079 {
1080  ASSERTL0(false, "This method is not defined for this expansion");
1081  StdExpansionSharedPtr returnval;
1082  return returnval;
1083 }
1084 
1086 {
1087  ASSERTL0(false, "This function has not been defined for this expansion");
1088  return false;
1089 }
1090 
1092 {
1093  return false;
1094 }
1095 
1097  const int dir, const Array<OneD, const NekDouble> &inarray,
1098  Array<OneD, NekDouble> &outarray)
1099 {
1100  boost::ignore_unused(dir, inarray, outarray);
1101  NEKERROR(ErrorUtil::efatal, "This method has not been defined");
1102 }
1103 
1104 /**
1105  *
1106  */
1108  const Array<OneD, const NekDouble> &direction,
1109  const Array<OneD, const NekDouble> &inarray,
1110  Array<OneD, NekDouble> &outarray)
1111 {
1112  boost::ignore_unused(direction, inarray, outarray);
1113  NEKERROR(ErrorUtil::efatal, "This method has not been defined");
1114 }
1115 
1116 /**
1117  *
1118  */
1120  const Array<OneD, const NekDouble> &inarray,
1121  Array<OneD, NekDouble> &outarray)
1122 {
1123  boost::ignore_unused(inarray, outarray);
1124  NEKERROR(ErrorUtil::efatal, "This method has not been defined");
1125 }
1126 
1127 /**
1128  * @brief Integrates the specified function over the domain.
1129  * @see StdRegions#StdExpansion#Integral.
1130  */
1132 {
1133  boost::ignore_unused(inarray);
1134  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1135  "local expansions");
1136  return 0;
1137 }
1138 
1139 /**
1140  * @brief Calculate the derivative of the physical points
1141  * @see StdRegions#StdExpansion#PhysDeriv
1142  */
1144  Array<OneD, NekDouble> &out_d1,
1145  Array<OneD, NekDouble> &out_d2,
1146  Array<OneD, NekDouble> &out_d3)
1147 {
1148  boost::ignore_unused(inarray, out_d1, out_d2, out_d3);
1149  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1150  "local expansions");
1151 }
1152 
1154  Array<OneD, NekDouble> &out_ds)
1155 {
1156  boost::ignore_unused(inarray, out_ds);
1157  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1158  "local expansions");
1159 }
1161  Array<OneD, NekDouble> &out_dn)
1162 {
1163  boost::ignore_unused(inarray, out_dn);
1164  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1165  "local expansions");
1166 }
1167 
1168 /**
1169  * @brief Calculate the derivative of the physical points in a
1170  * given direction
1171  * @see StdRegions#StdExpansion#PhysDeriv
1172  */
1173 void StdExpansion::v_PhysDeriv(const int dir,
1174  const Array<OneD, const NekDouble> &inarray,
1175  Array<OneD, NekDouble> &out_d0)
1176 
1177 {
1178  boost::ignore_unused(dir, inarray, out_d0);
1179  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1180  "specific element types");
1181 }
1182 
1183 /**
1184  * @brief Physical derivative along a direction vector.
1185  * @see StdRegions#StdExpansion#PhysDirectionalDeriv
1186  */
1188  const Array<OneD, const NekDouble> &inarray,
1189  const Array<OneD, const NekDouble> &direction,
1190  Array<OneD, NekDouble> &outarray)
1191 {
1192  boost::ignore_unused(inarray, direction, outarray);
1193  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1194  "specific element types");
1195 }
1196 
1198  Array<OneD, NekDouble> &out_d1,
1199  Array<OneD, NekDouble> &out_d2,
1200  Array<OneD, NekDouble> &out_d3)
1201 {
1202  boost::ignore_unused(inarray, out_d1, out_d2, out_d3);
1203  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1204 }
1205 
1206 void StdExpansion::v_StdPhysDeriv(const int dir,
1207  const Array<OneD, const NekDouble> &inarray,
1208  Array<OneD, NekDouble> &outarray)
1209 {
1210  boost::ignore_unused(dir, inarray, outarray);
1211  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1212 }
1213 
1215  const Array<OneD, const NekDouble> &coords,
1216  const Array<OneD, const NekDouble> &physvals)
1217 {
1218  boost::ignore_unused(coords, physvals);
1219  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1220  return 0;
1221 }
1222 
1225  const Array<OneD, const NekDouble> &physvals)
1226 {
1227  boost::ignore_unused(I, physvals);
1228  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1229  return 0;
1230 }
1231 
1233  const Array<OneD, const NekDouble> &coords, int mode)
1234 {
1235  boost::ignore_unused(coords, mode);
1236  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1237  return 0;
1238 }
1239 
1241  const Array<OneD, NekDouble> &coord,
1242  const Array<OneD, const NekDouble> &inarray,
1243  std::array<NekDouble, 3> &firstOrderDerivs)
1244 {
1245  boost::ignore_unused(coord, inarray, firstOrderDerivs);
1247  "PhysEvaluate first order derivative method does not exist"
1248  " for this shape type: " +
1249  static_cast<std::string>(
1251  return 0;
1252 }
1253 
1255  const Array<OneD, NekDouble> &coord,
1256  const Array<OneD, const NekDouble> &inarray,
1257  std::array<NekDouble, 3> &firstOrderDerivs,
1258  std::array<NekDouble, 6> &secondOrderDerivs)
1259 {
1260  boost::ignore_unused(coord, inarray, firstOrderDerivs, secondOrderDerivs);
1262  "PhysEvaluate second order derivative method does not exist"
1263  " for this shape type: " +
1264  static_cast<std::string>(
1266  return 0;
1267 }
1268 
1269 void StdExpansion::v_FillMode(const int mode, Array<OneD, NekDouble> &outarray)
1270 {
1271  boost::ignore_unused(mode, outarray);
1272  NEKERROR(ErrorUtil::efatal, "This function has not "
1273  "been defined for this shape");
1274 }
1275 
1277 {
1278  boost::ignore_unused(mkey);
1279  NEKERROR(ErrorUtil::efatal, "This function has not "
1280  "been defined for this element");
1281  DNekMatSharedPtr returnval;
1282  return returnval;
1283 }
1284 
1286 {
1287  boost::ignore_unused(mkey);
1288  NEKERROR(ErrorUtil::efatal, "This function has not "
1289  "been defined for this element");
1290  DNekMatSharedPtr returnval;
1291  return returnval;
1292 }
1293 
1295  Array<OneD, NekDouble> &coords_1,
1296  Array<OneD, NekDouble> &coords_2)
1297 {
1298  boost::ignore_unused(coords_0, coords_1, coords_2);
1299  NEKERROR(ErrorUtil::efatal, "Write coordinate definition method");
1300 }
1301 
1303  Array<OneD, NekDouble> &coord)
1304 {
1305  boost::ignore_unused(Lcoord, coord);
1306  NEKERROR(ErrorUtil::efatal, "Write coordinate definition method");
1307 }
1308 
1310 {
1311  return GetShapeDimension();
1312 }
1313 
1315 {
1316  boost::ignore_unused(outarray);
1317  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1318 }
1319 
1321 {
1322  boost::ignore_unused(outarray);
1323  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1324 }
1325 
1326 int StdExpansion::v_GetVertexMap(const int localVertexId, bool useCoeffPacking)
1327 {
1328  boost::ignore_unused(localVertexId, useCoeffPacking);
1329  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1330  return 0;
1331 }
1332 
1334  Array<OneD, unsigned int> &maparray,
1335  Array<OneD, int> &signarray,
1336  Orientation traceOrient, int P, int Q)
1337 {
1338  boost::ignore_unused(tid, maparray, signarray, traceOrient, P, Q);
1339  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1340 }
1341 
1342 void StdExpansion::v_GetTraceCoeffMap(const unsigned int traceid,
1343  Array<OneD, unsigned int> &maparray)
1344 {
1345  boost::ignore_unused(traceid, maparray);
1346  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1347 }
1348 
1349 void StdExpansion::v_GetElmtTraceToTraceMap(const unsigned int tid,
1350  Array<OneD, unsigned int> &maparray,
1351  Array<OneD, int> &signarray,
1352  Orientation traceOrient, int P,
1353  int Q)
1354 {
1355  boost::ignore_unused(tid, maparray, signarray, traceOrient, P, Q);
1356  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1357 }
1358 
1360  const int tid, Array<OneD, unsigned int> &maparray,
1361  Array<OneD, int> &signarray, const Orientation traceOrient)
1362 {
1363  boost::ignore_unused(tid, maparray, signarray, traceOrient);
1364  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1365 }
1366 
1367 void StdExpansion::v_GetTraceNumModes(const int tid, int &numModes0,
1368  int &numModes1, Orientation traceOrient)
1369 {
1370  boost::ignore_unused(tid, traceOrient, numModes0, numModes1);
1371  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1372 }
1373 
1375  const int vertex, const Array<OneD, const NekDouble> &inarray,
1376  NekDouble &outarray)
1377 {
1378  boost::ignore_unused(vertex, inarray, outarray);
1379  NEKERROR(ErrorUtil::efatal, "Method does not exist for "
1380  "this shape or library");
1381 }
1382 
1384  const Array<OneD, const NekDouble> &inarray,
1385  Array<OneD, NekDouble> &outarray)
1386 {
1387  v_MultiplyByStdQuadratureMetric(inarray, outarray);
1388 }
1389 
1391  const Array<OneD, const NekDouble> &inarray,
1392  Array<OneD, NekDouble> &outarray)
1393 {
1394  boost::ignore_unused(inarray, outarray);
1396  "Method does not exist for this shape or library");
1397 }
1398 
1400  const Array<OneD, const NekDouble> &inarray,
1401  Array<OneD, NekDouble> &outarray)
1402 {
1403  boost::ignore_unused(inarray, outarray);
1404  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1405 }
1406 
1408  const Array<OneD, const NekDouble> &inarray,
1409  Array<OneD, NekDouble> &outarray, bool multiplybyweights)
1410 {
1411  boost::ignore_unused(inarray, outarray, multiplybyweights);
1412  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1413 }
1414 
1415 /**
1416  *
1417  */
1419  const Array<OneD, const NekDouble> &direction,
1420  const Array<OneD, const NekDouble> &inarray,
1421  Array<OneD, NekDouble> &outarray)
1422 {
1423  boost::ignore_unused(direction, inarray, outarray);
1424  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1425 }
1426 
1428  const int dir, const Array<OneD, const NekDouble> &inarray,
1429  Array<OneD, NekDouble> &outarray)
1430 {
1431  boost::ignore_unused(dir, inarray, outarray);
1432  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1433 }
1434 
1436  Array<OneD, NekDouble> &outarray,
1437  const StdMatrixKey &mkey)
1438 {
1439  // If this function is not reimplemented on shape level, the function
1440  // below will be called
1441  MassMatrixOp_MatFree(inarray, outarray, mkey);
1442 }
1443 
1445  const Array<OneD, const NekDouble> &inarray,
1446  Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1447 {
1448  // If this function is not reimplemented on shape level, the function
1449  // below will be called
1450  LaplacianMatrixOp_MatFree(inarray, outarray, mkey);
1451 }
1452 
1454  const StdMatrixKey &mkey)
1455 {
1456  boost::ignore_unused(array, mkey);
1457  ASSERTL0(false, "This function is not defined in StdExpansion.");
1458 }
1459 
1461  const NekDouble alpha,
1462  const NekDouble exponent,
1463  const NekDouble cutoff)
1464 {
1465  boost::ignore_unused(array, alpha, exponent, cutoff);
1466  ASSERTL0(false, "This function is not defined in StdExpansion.");
1467 }
1468 
1470  int numMin, const Array<OneD, const NekDouble> &inarray,
1471  Array<OneD, NekDouble> &outarray)
1472 {
1473  boost::ignore_unused(numMin, inarray, outarray);
1474  ASSERTL0(false, "This function is not defined in StdExpansion.");
1475 }
1476 
1478  const int k1, const int k2, const Array<OneD, const NekDouble> &inarray,
1479  Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1480 {
1481  // If this function is not reimplemented on shape level, the function
1482  // below will be called
1483  LaplacianMatrixOp_MatFree(k1, k2, inarray, outarray, mkey);
1484 }
1485 
1487  const int i, const Array<OneD, const NekDouble> &inarray,
1488  Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1489 {
1490  // If this function is not reimplemented on shape level, the function
1491  // below will be called
1492  WeakDerivMatrixOp_MatFree(i, inarray, outarray, mkey);
1493 }
1494 
1496  const Array<OneD, const NekDouble> &inarray,
1497  Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1498 {
1499  // If this function is not reimplemented on shape level, the function
1500  // below will be called
1501  WeakDirectionalDerivMatrixOp_MatFree(inarray, outarray, mkey);
1502 }
1503 
1505  const Array<OneD, const NekDouble> &inarray,
1506  Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1507 {
1508  // If this function is not reimplemented on shape level, the function
1509  // below will be called
1510  MassLevelCurvatureMatrixOp_MatFree(inarray, outarray, mkey);
1511 }
1512 
1514  const Array<OneD, const NekDouble> &inarray,
1515  Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey,
1516  bool addDiffusionTerm)
1517 {
1518  // If this function is not reimplemented on shape level, the function
1519  // below will be called
1520  LinearAdvectionDiffusionReactionMatrixOp_MatFree(inarray, outarray, mkey,
1521  addDiffusionTerm);
1522 }
1523 
1525  const Array<OneD, const NekDouble> &inarray,
1526  Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1527 {
1528  // If this function is not reimplemented on shape level, the function
1529  // below will be called
1530  HelmholtzMatrixOp_MatFree(inarray, outarray, mkey);
1531 }
1532 
1534  const Array<OneD, const NekDouble> &inarray,
1535  Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1536 {
1537  // If this function is not reimplemented on shape level, the function
1538  // below will be called
1539  LaplacianMatrixOp_MatFree_GenericImpl(inarray, outarray, mkey);
1540 }
1541 
1543  const Array<OneD, const NekDouble> &inarray,
1545 {
1546  boost::ignore_unused(inarray, outarray, wsp);
1547  ASSERTL0(false, "Not implemented.");
1548 }
1549 
1551  const Array<OneD, const NekDouble> &inarray,
1552  Array<OneD, NekDouble> &outarray, const StdMatrixKey &mkey)
1553 {
1554  // If this function is not reimplemented on shape level, the function
1555  // below will be called
1556  HelmholtzMatrixOp_MatFree_GenericImpl(inarray, outarray, mkey);
1557 }
1558 
1560  const DNekScalMatSharedPtr &m_transformationmatrix)
1561 {
1562  boost::ignore_unused(m_transformationmatrix);
1563  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
1564  return NullDNekMatSharedPtr;
1565 }
1566 
1568  const Array<OneD, const NekDouble> &inarray,
1569  Array<OneD, NekDouble> &outarray, int npset)
1570 {
1572  DNekMatSharedPtr intmat;
1573 
1574  int nqtot = GetTotPoints();
1575  int np = 0;
1576  if (npset == -1) // use values from basis num points()
1577  {
1578  int nqbase;
1579  for (int i = 0; i < m_base.size(); ++i)
1580  {
1581  nqbase = m_base[i]->GetNumPoints();
1582  np = std::max(np, nqbase);
1583  }
1584 
1585  StdMatrixKey Ikey(ePhysInterpToEquiSpaced, shape, *this);
1586  intmat = GetStdMatrix(Ikey);
1587  }
1588  else
1589  {
1590  np = npset;
1591 
1592  ConstFactorMap cmap;
1593  cmap[eFactorConst] = np;
1594  StdMatrixKey Ikey(ePhysInterpToEquiSpaced, shape, *this, cmap);
1595  intmat = GetStdMatrix(Ikey);
1596  }
1597 
1598  NekVector<NekDouble> in(nqtot, inarray, eWrapper);
1600  LibUtilities::GetNumberOfCoefficients(shape, np, np, np), outarray,
1601  eWrapper);
1602  out = (*intmat) * in;
1603 }
1604 
1606  bool standard)
1607 {
1608  boost::ignore_unused(conn, standard);
1610  "GetSimplexEquiSpacedConnectivity not"
1611  " implemented for " +
1612  static_cast<std::string>(
1614 }
1615 
1617  const Array<OneD, const NekDouble> &inarray,
1618  Array<OneD, NekDouble> &outarray)
1619 {
1621 
1622  // inarray has to be consistent with NumModes definition
1623  // There is also a check in GetStdMatrix to see if all
1624  // modes are of the same size
1625  ConstFactorMap cmap;
1626 
1627  cmap[eFactorConst] = m_base[0]->GetNumModes();
1628  StdMatrixKey Ikey(eEquiSpacedToCoeffs, shape, *this, cmap);
1629  DNekMatSharedPtr intmat = GetStdMatrix(Ikey);
1630 
1631  NekVector<NekDouble> in(m_ncoeffs, inarray, eWrapper);
1632  NekVector<NekDouble> out(m_ncoeffs, outarray, eWrapper);
1633  out = (*intmat) * in;
1634 }
1635 
1636 } // namespace StdRegions
1637 } // 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
Defines a specification for a set of points.
Definition: Points.h:59
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
The base class for all shapes.
Definition: StdExpansion.h:71
void GetBoundaryMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:675
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual int v_GetVertexMap(int localVertexId, bool useCoeffPacking=false)
virtual void v_LocCollapsedToLocCoord(const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi)
virtual ~StdExpansion()
Destructor.
virtual void v_PhysDeriv_n(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn)
virtual void v_PhysDirectionalDeriv(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &outarray)
Physical derivative along a direction vector.
StdExpansion()
Default Constructor.
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:130
virtual void v_GetTraceToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140
virtual void v_GetCoord(const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
void GeneralMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void FillMode(const int mode, Array< OneD, NekDouble > &outarray)
This function fills the array outarray with the mode-th mode of the expansion.
Definition: StdExpansion.h:497
virtual NekDouble v_StdPhysEvaluate(const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
void PhysInterpToSimplexEquiSpaced(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int npset=-1)
This function performs an interpolation from the physical space points provided at input into an arra...
void WeakDirectionalDerivMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:807
void WeakDerivMatrixOp_MatFree(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
DNekBlkMatSharedPtr CreateStdStaticCondMatrix(const StdMatrixKey &mkey)
Create the static condensation of a matrix when using a boundary interior decomposition.
void LaplacianMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_IProductWRTDirectionalDerivBase_SumFac(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void MassLevelCurvatureMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:814
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_0, Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2)
void LinearAdvectionDiffusionReactionMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:609
virtual void v_GetElmtTraceToTraceMap(const unsigned int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
virtual void v_NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
virtual void v_GetSimplexEquiSpacedConnectivity(Array< OneD, int > &conn, bool standard=true)
virtual void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
virtual std::shared_ptr< StdExpansion > v_GetLinStdExp(void) const
void WeakDerivMatrixOp(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:799
virtual void v_GetInteriorMap(Array< OneD, unsigned int > &outarray)
virtual void v_SetCoeffsToOrientation(StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
NekDouble Linf(const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
Function to evaluate the discrete error where is given by the array sol.
void EquiSpacedToCoeffs(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs a projection/interpolation from the equispaced points sometimes used in post-p...
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:758
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_WeakDirectionalDerivMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual NekDouble v_PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
virtual DNekMatSharedPtr v_BuildInverseTransformationMatrix(const DNekScalMatSharedPtr &m_transformationmatrix)
NekDouble L2(const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
Function to evaluate the discrete error, where is given by the array sol.
virtual void v_GetVertexPhysVals(const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray)
void IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
virtual void v_MassLevelCurvatureMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_GetTraceCoeffMap(const unsigned int traceid, Array< OneD, unsigned int > &maparray)
virtual bool v_IsBoundaryInteriorExpansion() const
virtual void v_MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:765
virtual void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray)
void LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_ExponentialFilter(Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
void MassLevelCurvatureMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:373
virtual void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey)
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:680
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey)
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_GetTraceNumModes(const int fid, int &numModes0, int &numModes1, Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2)
virtual void v_PhysDeriv_s(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds)
virtual void v_FwdTransBndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray)
void GenStdMatBwdDeriv(const int dir, DNekMatSharedPtr &mat)
virtual void v_MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual const LibUtilities::PointsKey v_GetNodalPointsKey() const
virtual void v_WeakDerivMatrixOp(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LinearAdvectionDiffusionReactionMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
Definition: StdExpansion.h:821
void IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
virtual void v_LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
virtual NekDouble v_PhysEvaluateBasis(const Array< OneD, const NekDouble > &coords, int mode)
virtual void v_GetTraceInteriorToElementMap(const int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eForwards)
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)=0
Calculates the inner product of a given function f with the different modes of the expansion.
virtual void v_LinearAdvectionDiffusionReactionMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
void HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:836
virtual void v_ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_DropLocStaticCondMatrix(const LocalRegions::MatrixKey &mkey)
void WeakDirectionalDerivMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_LaplacianMatrixOp_MatFree_Kernel(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
NekDouble H1(const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
Function to evaluate the discrete error, where is given by the array sol.
virtual int v_CalcNumberOfCoefficients(const std::vector< unsigned int > &nummodes, int &modes_offset)
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray)
Integrates the specified function over the domain.
void SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
Definition: StdExpansion.h:779
virtual std::shared_ptr< StdExpansion > v_GetStdExp() const
void HelmholtzMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)=0
virtual LibUtilities::PointsKey v_GetTracePointsKey(const int i, const int j) const
Array< OneD, LibUtilities::BasisSharedPtr > m_base
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
virtual void v_IProductWRTDirectionalDerivBase(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2, Array< OneD, NekDouble > &out_d3)
Calculate the derivative of the physical points.
virtual const LibUtilities::BasisKey v_GetTraceBasisKey(const int i, const int k) const
virtual void v_StdPhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2, Array< OneD, NekDouble > &out_d3)
NekDouble StdPhysEvaluate(const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
void MassMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void GeneralMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix(const LocalRegions::MatrixKey &mkey)
const ConstFactorMap & GetConstFactors() const
Definition: StdMatrixKey.h:140
LibUtilities::ShapeType GetShapeType() const
Definition: StdMatrixKey.h:90
const VarCoeffMap & GetVarCoeffs() const
Definition: StdMatrixKey.h:171
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:85
const Array< OneD, const NekDouble > & GetVarCoeff(const StdRegions::VarCoeffType &coeff) const
Definition: StdMatrixKey.h:150
bool HasVarCoeff(const StdRegions::VarCoeffType &coeff) const
Definition: StdMatrixKey.h:176
LibUtilities::PointsType GetNodalPointsType() const
Definition: StdMatrixKey.h:95
NekDouble GetConstFactor(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:126
bool ConstFactorExists(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:135
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
Definition: Blas.hpp:154
const char *const ShapeTypeMap[SIZE_ShapeType]
Definition: ShapeType.hpp:79
BasisManagerT & BasisManager(void)
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.
int GetNumberOfCoefficients(ShapeType shape, std::vector< unsigned int > &modes, int offset=0)
Definition: ShapeType.hpp:310
static const PointsKey NullPointsKey(0, eNoPointsType)
static const NekDouble kNekSqrtTol
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
const char *const MatrixTypeMap[]
Definition: StdRegions.hpp:143
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:399
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:400
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:344
std::shared_ptr< StdMatrixKey > StdMatrixKeySharedPtr
Definition: StdMatrixKey.h:201
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:77
static DNekScalBlkMatSharedPtr NullDNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:85
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:79
static DNekMatSharedPtr NullDNekMatSharedPtr
Definition: NekTypeDefs.hpp:83
static Array< OneD, NekDouble > NullNekDouble1DArray
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
void Vabs(int n, const T *x, const int incx, T *y, const int incy)
vabs: y = |x|
Definition: Vmath.cpp:553
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 Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:359
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:248
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492
T Vamax(int n, const T *x, const int incx)
Return the maximum absolute element in x called vamax to avoid conflict with max.
Definition: Vmath.cpp:999
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:419
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294