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 
39 #include <LibUtilities/Foundations/ManagerAccess.h> // for BasisManager, etc
40 
41 namespace Nektar
42 {
43  namespace StdRegions
44  {
46  m_elmt_id(0),
47  m_ncoeffs(0)
48  {
49  }
50 
51  StdExpansion::StdExpansion(const int numcoeffs, const int numbases,
52  const LibUtilities::BasisKey &Ba,
53  const LibUtilities::BasisKey &Bb,
54  const LibUtilities::BasisKey &Bc):
55  m_base(numbases),
56  m_elmt_id(0),
57  m_ncoeffs(numcoeffs),
59  std::bind(&StdExpansion::CreateStdMatrix, this, std::placeholders::_1),
60  std::string("StdExpansionStdMatrix")),
62  std::bind(&StdExpansion::CreateStdStaticCondMatrix, this, std::placeholders::_1),
63  std::string("StdExpansionStdStaticCondMatrix")),
65  std::bind(&StdExpansion::CreateIndexMap,this, std::placeholders::_1),
66  std::string("StdExpansionIndexMap"))
67  {
68  switch(m_base.num_elements())
69  {
70  case 3:
72  "NULL Basis attempting to be used.");
74  /* Falls through. */
75  case 2:
77  "NULL Basis attempting to be used.");
79  /* Falls through. */
80  case 1:
82  "NULL Basis attempting to be used.");
84  break;
85  default:
86  break;
87 // ASSERTL0(false, "numbases incorrectly specified");
88  };
89 
90  } //end constructor
91 
92 
94  std::enable_shared_from_this<StdExpansion>(T),
95  m_base(T.m_base),
101  {
102  }
103 
105  {
106  }
107 
109  const Array<OneD, const NekDouble>& sol)
110  {
111  NekDouble val;
112  int ntot = GetTotPoints();
113  Array<OneD, NekDouble> wsp(ntot);
114 
115  if(sol == NullNekDouble1DArray)
116  {
117  Vmath::Vabs(ntot, phys, 1, wsp, 1);
118  }
119  else
120  {
121  Vmath::Vsub(ntot, sol, 1, phys, 1, wsp, 1);
122  Vmath::Vabs(ntot, wsp, 1, wsp, 1);
123  }
124 
125  val = Vmath::Vamax(ntot, wsp, 1);
126 
127  return val;
128  }
129 
131  const Array<OneD, const NekDouble>& sol)
132  {
133  NekDouble val;
134  int ntot = GetTotPoints();
135  Array<OneD, NekDouble> wsp(ntot);
136 
137  if (sol.num_elements() == 0)
138  {
139  Vmath::Vmul(ntot, phys, 1, phys, 1, wsp, 1);
140  }
141  else
142  {
143  Vmath::Vsub(ntot, sol, 1, phys, 1, wsp, 1);
144  Vmath::Vmul(ntot, wsp, 1, wsp, 1, wsp, 1);
145  }
146 
147  val = v_Integral(wsp);
148 
149  // if val too small, sqrt returns nan.
151  {
152  return 0.0;
153  }
154  else
155  {
156  return sqrt(val);
157  }
158  }
159 
161  const Array<OneD, const NekDouble>& sol)
162  {
163  int i;
164  NekDouble val;
165  int ntot = GetTotPoints();
166  int coordim = v_GetCoordim();
167  Array<OneD, NekDouble> wsp(3*ntot);
168  Array<OneD, NekDouble> wsp_deriv = wsp + ntot;
169  Array<OneD, NekDouble> sum = wsp_deriv + ntot;
170 
171  if(sol == NullNekDouble1DArray)
172  {
173  Vmath::Vcopy(ntot,phys, 1, wsp, 1);
174  Vmath::Vmul(ntot, phys, 1, phys, 1, sum, 1);
175  }
176  else
177  {
178  Vmath::Vsub(ntot, sol, 1, phys, 1, wsp, 1);
179  Vmath::Vmul(ntot, wsp, 1, wsp, 1, sum, 1);
180  }
181 
182 
183  for(i = 0; i < coordim; ++i)
184  {
185  v_PhysDeriv(i,wsp,wsp_deriv);
186  Vmath::Vvtvp(ntot,wsp_deriv,1,wsp_deriv,1,sum,1,sum,1);
187  }
188 
189  val = sqrt(v_Integral(sum));
190 
191  return val;
192  }
193 
194 
196  {
197  DNekBlkMatSharedPtr returnval;
198 
199  DNekMatSharedPtr mat = GetStdMatrix(mkey);
200  int nbdry = NumBndryCoeffs(); // also checks to see if this is a boundary interior decomposed expansion
201  int nint = m_ncoeffs - nbdry;
206 
207  int i,j;
208 
209  Array<OneD,unsigned int> bmap(nbdry);
210  Array<OneD,unsigned int> imap(nint);
211  GetBoundaryMap(bmap);
212  GetInteriorMap(imap);
213 
214  for(i = 0; i < nbdry; ++i)
215  {
216  for(j = 0; j < nbdry; ++j)
217  {
218  (*A)(i,j) = (*mat)(bmap[i],bmap[j]);
219  }
220 
221  for(j = 0; j < nint; ++j)
222  {
223  (*B)(i,j) = (*mat)(bmap[i],imap[j]);
224  }
225  }
226 
227  for(i = 0; i < nint; ++i)
228  {
229  for(j = 0; j < nbdry; ++j)
230  {
231  (*C)(i,j) = (*mat)(imap[i],bmap[j]);
232  }
233 
234  for(j = 0; j < nint; ++j)
235  {
236  (*D)(i,j) = (*mat)(imap[i],imap[j]);
237  }
238  }
239 
240  // Calculate static condensed system
241  if(nint)
242  {
243  D->Invert();
244  (*B) = (*B)*(*D);
245  (*A) = (*A) - (*B)*(*C);
246  }
247 
248  // set up block matrix system
249  Array<OneD, unsigned int> exp_size(2);
250  exp_size[0] = nbdry;
251  exp_size[1] = nint;
252  returnval = MemoryManager<DNekBlkMat>::AllocateSharedPtr(exp_size,exp_size);
253 
254  returnval->SetBlock(0,0,A);
255  returnval->SetBlock(0,1,B);
256  returnval->SetBlock(1,0,C);
257  returnval->SetBlock(1,1,D);
258 
259  return returnval;
260  }
261 
263  {
264  IndexMapValuesSharedPtr returnval;
265 
266  IndexMapType itype = ikey.GetIndexMapType();
267 
268  int entity = ikey.GetIndexEntity();
269 
270  Orientation orient = ikey.GetIndexOrientation();
271 
274 
275  switch(itype)
276  {
277  case eEdgeToElement:
278  {
279  v_GetEdgeToElementMap(entity,orient,map,sign);
280  }
281  break;
282  case eFaceToElement:
283  {
284  v_GetFaceToElementMap(entity,orient,map,sign);
285  }
286  break;
287  case eEdgeInterior:
288  {
289  v_GetEdgeInteriorMap(entity,orient,map,sign);
290  }
291  break;
292  case eFaceInterior:
293  {
294  v_GetFaceInteriorMap(entity,orient,map,sign);
295  }
296  break;
297  case eBoundary:
298  {
299  ASSERTL0(false,"Boundary Index Map not implemented yet.");
300  }
301  break;
302  case eVertex:
303  {
304  ASSERTL0(false,"Vertex Index Map not implemented yet.");
305  }
306  break;
307  default:
308  {
309  ASSERTL0(false,"The Index Map you are requiring is not between the possible options.");
310  }
311  }
312 
313  returnval = MemoryManager<IndexMapValues>::AllocateSharedPtr(map.num_elements());
314 
315  for(int i = 0; i < map.num_elements(); i++)
316  {
317  (*returnval)[i].index = map[i];
318  (*returnval)[i].sign = sign[i];
319  }
320 
321  return returnval;
322  }
323 
325  {
326  int i;
327  DNekMatSharedPtr returnval;
328 
329  switch(mkey.GetMatrixType())
330  {
331  case eInvMass:
332  {
334  DNekMatSharedPtr mmat = GetStdMatrix(masskey);
335 
336  returnval = MemoryManager<DNekMat>::AllocateSharedPtr(*mmat); //Populate standard mass matrix.
337  returnval->Invert();
338  }
339  break;
340  case eInvNBasisTrans:
341  {
343  DNekMatSharedPtr tmpmat = GetStdMatrix(tmpkey);
344  returnval = MemoryManager<DNekMat>::AllocateSharedPtr(*tmpmat); //Populate matrix.
345  returnval->Invert();
346  }
347  break;
348  case eBwdTrans:
349  {
350  int nq = GetTotPoints();
352  Array<OneD, NekDouble> tmpout(nq);
353 
355 
356  for(int i=0; i<m_ncoeffs; ++i)
357  {
358  Vmath::Zero(m_ncoeffs, tmpin, 1);
359  tmpin[i] = 1.0;
360 
361  BwdTrans_SumFac(tmpin,tmpout);
362 
363  Vmath::Vcopy(nq,tmpout.get(),1,
364  returnval->GetRawPtr()+i*nq,1);
365  }
366  }
367  break;
368  case eIProductWRTBase:
369  {
370  int nq = GetTotPoints();
371  Array<OneD, NekDouble> tmpin(nq);
373 
375 
376  for(i=0; i < nq; ++i)
377  {
378  Vmath::Zero(nq, tmpin, 1);
379  tmpin[i] = 1.0;
380 
381  IProductWRTBase_SumFac(tmpin,tmpout);
382 
383  Vmath::Vcopy(m_ncoeffs,tmpout.get(),1,
384  returnval->GetRawPtr()+i*m_ncoeffs,1);
385  }
386  }
387  break;
389  {
390  int nq = GetTotPoints();
391  Array<OneD, NekDouble> tmpin(nq);
393 
395 
396  for(i=0; i < nq; ++i)
397  {
398  Vmath::Zero(nq, tmpin, 1);
399  tmpin[i] = 1.0;
400 
401  IProductWRTDerivBase_SumFac(0,tmpin,tmpout);
402 
403  Vmath::Vcopy(m_ncoeffs,tmpout.get(),1,
404  returnval->GetRawPtr()+i*m_ncoeffs,1);
405  }
406  }
407  break;
409  {
410  int nq = GetTotPoints();
411  Array<OneD, NekDouble> tmpin(nq);
413 
415 
416  for(i=0; i < nq; ++i)
417  {
418  Vmath::Zero(nq, tmpin, 1);
419  tmpin[i] = 1.0;
420 
421  IProductWRTDerivBase_SumFac(1,tmpin,tmpout);
422 
423  Vmath::Vcopy(m_ncoeffs,tmpout.get(),1,
424  returnval->GetRawPtr()+i*m_ncoeffs,1);
425  }
426  }
427  break;
429  {
430  int nq = GetTotPoints();
431  Array<OneD, NekDouble> tmpin(nq);
433 
435 
436  for(i=0; i < nq; ++i)
437  {
438  Vmath::Zero(nq, tmpin, 1);
439  tmpin[i] = 1.0;
440 
441  IProductWRTDerivBase_SumFac(2,tmpin,tmpout);
442 
443  Vmath::Vcopy(m_ncoeffs,tmpout.get(),1,
444  returnval->GetRawPtr()+i*m_ncoeffs,1);
445  }
446  }
447  break;
448  case eEquiSpacedToCoeffs:
449  {
450  // check to see if equispaced basis
451  int nummodes = m_base[0]->GetNumModes();
452  bool equispaced = true;
453  for(int i = 1; i < m_base.num_elements(); ++i)
454  {
455  if(m_base[i]->GetNumModes() != nummodes)
456  {
457  equispaced = false;
458  }
459  }
460 
461  ASSERTL0(equispaced,
462  "Currently need to have same num modes in all "
463  "directionmodes to use EquiSpacedToCoeff method");
464 
465  int ntot = GetTotPoints();
466  Array<OneD, NekDouble> qmode(ntot);
468 
471  for(int i = 0; i < m_ncoeffs; ++i)
472  {
473  // Get mode at quadrature points
474  FillMode(i,qmode);
475 
476  // interpolate to equi spaced
477  PhysInterpToSimplexEquiSpaced(qmode,emode,nummodes);
478 
479  // fill matrix
480  Vmath::Vcopy(m_ncoeffs, &emode[0], 1,
481  returnval->GetRawPtr() + i*m_ncoeffs, 1);
482  }
483  // invert matrix
484  returnval->Invert();
485 
486  }
487  break;
488  case eMass:
489  case eHelmholtz:
490  case eLaplacian:
491  case eLaplacian00:
492  case eLaplacian01:
493  case eLaplacian02:
494  case eLaplacian11:
495  case eLaplacian12:
496  case eLaplacian22:
497  case eWeakDeriv0:
498  case eWeakDeriv1:
499  case eWeakDeriv2:
501  case eMassLevelCurvature:
504  {
507  DNekMat &Mat = *returnval;
508 
509  for(i=0; i < m_ncoeffs; ++i)
510  {
511  Vmath::Zero(m_ncoeffs, tmp, 1);
512  tmp[i] = 1.0;
513 
514  GeneralMatrixOp_MatFree(tmp,tmp,mkey);
515 
516  Vmath::Vcopy(m_ncoeffs,&tmp[0],1,
517  &(Mat.GetPtr())[0]+i*m_ncoeffs,1);
518  }
519  }
520  break;
521  default:
522  {
523  NEKERROR(ErrorUtil::efatal, "This type of matrix can not be created using a general approach");
524  }
525  break;
526  }
527 
528  return returnval;
529  }
530 
532  Array<OneD,NekDouble> &outarray,
533  const StdMatrixKey &mkey)
534  {
535  switch(mkey.GetMatrixType())
536  {
537  case eMass:
538  MassMatrixOp(inarray,outarray,mkey);
539  break;
540  case eWeakDeriv0:
541  WeakDerivMatrixOp(0,inarray,outarray,mkey);
542  break;
543  case eWeakDeriv1:
544  WeakDerivMatrixOp(1,inarray,outarray,mkey);
545  break;
546  case eWeakDeriv2:
547  WeakDerivMatrixOp(2,inarray,outarray,mkey);
548  break;
550  WeakDirectionalDerivMatrixOp(inarray,outarray,mkey);
551  break;
552  case eMassLevelCurvature:
553  MassLevelCurvatureMatrixOp(inarray,outarray,mkey);
554  break;
556  LinearAdvectionDiffusionReactionMatrixOp(inarray,outarray,mkey,false);
557  break;
559  LinearAdvectionDiffusionReactionMatrixOp(inarray,outarray,mkey);
560  break;
561  case eLaplacian:
562  LaplacianMatrixOp(inarray,outarray,mkey);
563  break;
564  case eLaplacian00:
565  LaplacianMatrixOp(0,0,inarray,outarray,mkey);
566  break;
567  case eLaplacian01:
568  LaplacianMatrixOp(0,1,inarray,outarray,mkey);
569  break;
570  case eLaplacian02:
571  LaplacianMatrixOp(0,2,inarray,outarray,mkey);
572  break;
573  case eLaplacian10:
574  LaplacianMatrixOp(1,0,inarray,outarray,mkey);
575  break;
576  case eLaplacian11:
577  LaplacianMatrixOp(1,1,inarray,outarray,mkey);
578  break;
579  case eLaplacian12:
580  LaplacianMatrixOp(1,2,inarray,outarray,mkey);
581  break;
582  case eLaplacian20:
583  LaplacianMatrixOp(2,0,inarray,outarray,mkey);
584  break;
585  case eLaplacian21:
586  LaplacianMatrixOp(2,1,inarray,outarray,mkey);
587  break;
588  case eLaplacian22:
589  LaplacianMatrixOp(2,2,inarray,outarray,mkey);
590  break;
591  case eHelmholtz:
592  HelmholtzMatrixOp(inarray,outarray,mkey);
593  break;
594  default:
595  NEKERROR(ErrorUtil::efatal, "This matrix does not have an operator");
596  break;
597  }
598  }
599 
601  Array<OneD,NekDouble> &outarray,
602  const StdMatrixKey &mkey)
603  {
604  switch(mkey.GetMatrixType())
605  {
606  case eMass:
607  MassMatrixOp_MatFree(inarray,outarray,mkey);
608  break;
609  case eWeakDeriv0:
610  WeakDerivMatrixOp_MatFree(0,inarray,outarray,mkey);
611  break;
612  case eWeakDeriv1:
613  WeakDerivMatrixOp_MatFree(1,inarray,outarray,mkey);
614  break;
615  case eWeakDeriv2:
616  WeakDerivMatrixOp_MatFree(2,inarray,outarray,mkey);
617  break;
619  WeakDirectionalDerivMatrixOp_MatFree(inarray,outarray,mkey);
620  break;
621  case eMassLevelCurvature:
622  MassLevelCurvatureMatrixOp_MatFree(inarray,outarray,mkey);
623  break;
625  LinearAdvectionDiffusionReactionMatrixOp_MatFree(inarray,outarray,mkey,false);
626  break;
629  break;
630  case eLaplacian:
631  LaplacianMatrixOp_MatFree(inarray,outarray,mkey);
632  break;
633  case eLaplacian00:
634  LaplacianMatrixOp_MatFree(0,0,inarray,outarray,mkey);
635  break;
636  case eLaplacian01:
637  LaplacianMatrixOp_MatFree(0,1,inarray,outarray,mkey);
638  break;
639  case eLaplacian02:
640  LaplacianMatrixOp_MatFree(0,2,inarray,outarray,mkey);
641  break;
642  case eLaplacian10:
643  LaplacianMatrixOp_MatFree(1,0,inarray,outarray,mkey);
644  break;
645  case eLaplacian11:
646  LaplacianMatrixOp_MatFree(1,1,inarray,outarray,mkey);
647  break;
648  case eLaplacian12:
649  LaplacianMatrixOp_MatFree(1,2,inarray,outarray,mkey);
650  break;
651  case eLaplacian20:
652  LaplacianMatrixOp_MatFree(2,0,inarray,outarray,mkey);
653  break;
654  case eLaplacian21:
655  LaplacianMatrixOp_MatFree(2,1,inarray,outarray,mkey);
656  break;
657  case eLaplacian22:
658  LaplacianMatrixOp_MatFree(2,2,inarray,outarray,mkey);
659  break;
660  case eHelmholtz:
661  HelmholtzMatrixOp_MatFree(inarray,outarray,mkey);
662  break;
663  default:
664  NEKERROR(ErrorUtil::efatal, "This matrix does not have an operator");
665  break;
666  }
667  }
668 
670  Array<OneD,NekDouble> &outarray,
671  const StdMatrixKey &mkey)
672  {
673  int nq = GetTotPoints();
674  Array<OneD, NekDouble> tmp(nq);
675 
676  v_BwdTrans(inarray,tmp);
677 
678  if(mkey.HasVarCoeff(eVarCoeffMass))
679  {
680  Vmath::Vmul(nq, mkey.GetVarCoeff(eVarCoeffMass), 1, tmp, 1, tmp, 1);
681  }
682 
683  v_IProductWRTBase(tmp, outarray);
684  }
685 
686  void StdExpansion::LaplacianMatrixOp_MatFree(const int k1, const int k2,
687  const Array<OneD, const NekDouble> &inarray,
688  Array<OneD,NekDouble> &outarray,
689  const StdMatrixKey &mkey)
690  {
691  ASSERTL1(k1 >= 0 && k1 < GetCoordim(),"invalid first argument");
692  ASSERTL1(k2 >= 0 && k2 < GetCoordim(),"invalid second argument");
693 
694  int nq = GetTotPoints();
695  Array<OneD, NekDouble> tmp(nq);
696  Array<OneD, NekDouble> dtmp(nq);
697  VarCoeffType varcoefftypes[3][3]
701  };
702 
703  v_BwdTrans(inarray,tmp);
704  v_PhysDeriv(k2,tmp,dtmp);
705  if (mkey.GetNVarCoeff()&&
707  {
708  if (k1 == k2)
709  {
710  // By default, k1 == k2 has \sigma = 1 (diagonal entries)
711  if(mkey.HasVarCoeff(varcoefftypes[k1][k1]))
712  {
713  Vmath::Vmul(nq, mkey.GetVarCoeff(varcoefftypes[k1][k1]), 1, dtmp, 1, dtmp, 1);
714  }
715  v_IProductWRTDerivBase(k1, dtmp, outarray);
716  }
717  else
718  {
719  // By default, k1 != k2 has \sigma = 0 (off-diagonal entries)
720  if(mkey.HasVarCoeff(varcoefftypes[k1][k2]))
721  {
722  Vmath::Vmul(nq, mkey.GetVarCoeff(varcoefftypes[k1][k2]), 1, dtmp, 1, dtmp, 1);
723  v_IProductWRTDerivBase(k1, dtmp, outarray);
724  }
725  else
726  {
727  Vmath::Zero(GetNcoeffs(), outarray, 1);
728  }
729  }
730  }
731  else
732  {
733  // Multiply by svv tensor
735  {
736  Vmath::Vcopy(nq, dtmp, 1, tmp, 1);
737  SVVLaplacianFilter(dtmp,mkey);
738  Vmath::Vadd(nq, tmp, 1, dtmp, 1, dtmp, 1);
739  }
740  v_IProductWRTDerivBase(k1, dtmp, outarray);
741  }
742  }
743 
745  Array<OneD,NekDouble> &outarray,
746  const StdMatrixKey &mkey)
747  {
748  const int dim = GetCoordim();
749 
750  int i,j;
751 
753  Array<OneD,NekDouble> store2(m_ncoeffs,0.0);
754 
755  if(mkey.GetNVarCoeff() == 0||mkey.ConstFactorExists(eFactorSVVDiffCoeff))
756  {
757  // just call diagonal matrix form of laplcian operator
758  for(i = 0; i < dim; ++i)
759  {
760  LaplacianMatrixOp(i,i,inarray,store,mkey);
761  Vmath::Vadd(m_ncoeffs, store, 1, store2, 1, store2, 1);
762  }
763  }
764  else
765  {
766  const MatrixType mtype[3][3]
770  StdMatrixKeySharedPtr mkeyij;
771 
772  for(i = 0; i < dim; i++)
773  {
774  for(j = 0; j < dim; j++)
775  {
776  mkeyij = MemoryManager<StdMatrixKey>::AllocateSharedPtr(mkey,mtype[i][j]);
777  LaplacianMatrixOp(i,j,inarray,store,*mkeyij);
778  Vmath::Vadd(m_ncoeffs, store, 1, store2, 1, store2, 1);
779  }
780  }
781  }
782 
783  Vmath::Vcopy(m_ncoeffs,store2.get(),1,outarray.get(),1);
784  }
785 
787  const Array<OneD, const NekDouble> &inarray,
788  Array<OneD,NekDouble> &outarray,
789  const StdMatrixKey &mkey)
790  {
792  int nq = GetTotPoints();
793 
794  v_BwdTrans(inarray,tmp);
795  v_PhysDeriv(k1,tmp,tmp);
796 
798  if(mkey.HasVarCoeff(keys[k1]))
799  {
800  Vmath::Vmul(nq, &(mkey.GetVarCoeff(keys[k1]))[0], 1, &tmp[0], 1, &tmp[0], 1);
801  }
802 
803  v_IProductWRTBase(tmp, outarray);
804  }
805 
807  const Array<OneD, const NekDouble> &inarray,
808  Array<OneD,NekDouble> &outarray,
809  const StdMatrixKey &mkey)
810  {
811  int nq = GetTotPoints();
812 
813  Array<OneD, NekDouble> tmp(nq), Dtmp(nq);
814  Array<OneD, NekDouble> Mtmp(nq), Mout(m_ncoeffs);
815 
816  v_BwdTrans(inarray,tmp);
818 
819  v_IProductWRTBase(Dtmp, outarray);
820 
821  // Compte M_{div tv}
822  Vmath::Vmul(nq, &(mkey.GetVarCoeff(eVarCoeffMFDiv))[0], 1,
823  &tmp[0], 1,
824  &Mtmp[0], 1);
825 
826  v_IProductWRTBase(Mtmp, Mout);
827 
828  // Add D_tv + M_{div tv}
829  Vmath::Vadd(m_ncoeffs, &Mout[0], 1,
830  &outarray[0], 1,
831  &outarray[0], 1);
832  }
833 
835  Array<OneD,NekDouble> &outarray,
836  const StdMatrixKey &mkey)
837  {
838  boost::ignore_unused(inarray, outarray, mkey);
839  ///@todo fix this
840  // int nqtot = GetTotPoints();
841  // int matrixid = mkey.GetMatrixID();
842  //
843  // NekDouble checkweight=0.0;
844  // Array<OneD, NekDouble> tmp(nqtot), tan(nqtot), dtan0(nqtot), dtan1(nqtot), weight(nqtot,0.0);
845  //
846  // int gmatnumber = (mkey.GetVariableCoefficient(1)).num_elements();
847  //
848  // v_BwdTrans(inarray,tmp);
849  //
850  // // weight = \grad \cdot tanvec
851  // for(int k = 0; k < GetCoordim(); ++k)
852  // {
853  // Vmath::Vcopy(nqtot, &(mkey.GetVariableCoefficient(0))[k*nqtot],
854  // 1, &tan[0], 1);
855  //
856  // // For Regular mesh ...
857  // if(gmatnumber==1)
858  // {
859  // // D_{/xi} and D_{/eta}
860  // v_PhysDeriv(0,tan,dtan0);
861  // v_PhysDeriv(1,tan,dtan1);
862  //
863  // // d v / d x_i = (d \xi / d x_i)*( d v / d \xi ) + (d \eta / d x_i)*( d v / d \eta )
864  // Vmath::Svtvp(nqtot,(mkey.GetVariableCoefficient(2*k+1))[0],&dtan0[0],1,&weight[0],1,&weight[0],1);
865  // Vmath::Svtvp(nqtot,(mkey.GetVariableCoefficient(2*k+2))[0],&dtan1[0],1,&weight[0],1,&weight[0],1);
866  // }
867  //
868  // // For Curved mesh ...
869  // else if(gmatnumber==nqtot)
870  // {
871  // // D_{x} and D_{y}
872  // v_PhysDeriv(k,tan,dtan0);
873  // Vmath::Vadd(nqtot,&dtan0[0],1,&weight[0],1,&weight[0],1);
874  // }
875  //
876  // else
877  // {
878  // ASSERTL1( ((gmatnumber=1) || (gmatnumber==nqtot) ), "Gmat is not in a right size");
879  // }
880  // }
881  //
882  // Vmath::Vmul(nqtot, &weight[0], 1, &tmp[0], 1, &tmp[0], 1);
883  // v_IProductWRTBase(tmp, outarray);
884  }
885 
887  Array<OneD,NekDouble> &outarray,
888  const StdMatrixKey &mkey,
889  bool addDiffusionTerm)
890  {
891 
892  int i;
893  int ndir = mkey.GetNVarCoeff(); // assume num.r consts corresponds to directions
894  ASSERTL0(ndir,"Must define at least one advection velocity");
895 
896  NekDouble lambda = mkey.GetConstFactor(eFactorLambda);
897  int totpts = GetTotPoints();
898  Array<OneD, NekDouble> tmp(3*totpts);
899  Array<OneD, NekDouble> tmp_deriv = tmp + totpts;
900  Array<OneD, NekDouble> tmp_adv = tmp_deriv + totpts;
901 
902 
903  ASSERTL1(ndir <= GetCoordim(),"Number of constants is larger than coordinate dimensions");
904 
905  v_BwdTrans(inarray,tmp);
906 
908 
909  //calculate u dx + v dy + ..
910  Vmath::Zero(totpts,tmp_adv,1);
911  for(i = 0; i < ndir; ++i)
912  {
913  v_PhysDeriv(i,tmp,tmp_deriv);
914  Vmath::Vvtvp(totpts,mkey.GetVarCoeff(varcoefftypes[i]),1,tmp_deriv,1,tmp_adv,1,tmp_adv,1);
915  }
916 
917  if(lambda) // add -lambda*u
918  {
919  Vmath::Svtvp(totpts,-lambda,tmp,1,tmp_adv,1,tmp_adv,1);
920  }
921 
922 
923  if(addDiffusionTerm)
924  {
926  StdMatrixKey mkeylap(eLaplacian,DetShapeType(),*this);
927  LaplacianMatrixOp(inarray,lap,mkeylap);
928 
929  v_IProductWRTBase(tmp_adv, outarray);
930  // Lap v - u.grad v + lambda*u
931  // => (grad u, grad v) + u.grad v - lambda*u
932  Vmath::Vadd(m_ncoeffs,lap,1,outarray,1,outarray,1);
933  }
934  else
935  {
936  v_IProductWRTBase(tmp_adv, outarray);
937  }
938 
939  }
940 
941 
943  Array<OneD,NekDouble> &outarray,
944  const StdMatrixKey &mkey)
945  {
946  NekDouble lambda = mkey.GetConstFactor(eFactorLambda);
948  StdMatrixKey mkeymass(eMass,DetShapeType(),*this);
949  StdMatrixKey mkeylap(eLaplacian,DetShapeType(),*this);
950 
951  MassMatrixOp(inarray,tmp,mkeymass);
952  LaplacianMatrixOp(inarray,outarray,mkeylap);
953 
954  Blas::Daxpy(m_ncoeffs, lambda, tmp, 1, outarray, 1);
955  }
956 
958  Array<OneD, NekDouble> &outarray)
959  {
960  int nq = GetTotPoints();
961  StdMatrixKey bwdtransmatkey(eBwdTrans,DetShapeType(),*this);
962  DNekMatSharedPtr bwdtransmat = GetStdMatrix(bwdtransmatkey);
963 
964  Blas::Dgemv('N',nq,m_ncoeffs,1.0,bwdtransmat->GetPtr().get(),
965  nq, inarray.get(), 1, 0.0, outarray.get(), 1);
966  }
967 
968  // VIRTUAL INLINE FUNCTIONS FROM HEADER FILE
969  void StdExpansion::SetUpPhysNormals(const int edge)
970  {
971  v_SetUpPhysNormals(edge);
972  }
973 
975  const Array<OneD, const NekDouble> &physvals)
976  {
977  return v_StdPhysEvaluate(Lcoord,physvals);
978  }
979 
981  {
982  return m_elmt_id;
983  }
984 
986  {
987  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
988  return NullNekDouble1DArray;
989  }
990 
991 
993  {
994  boost::ignore_unused(normal);
995  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
996  }
997 
998  void StdExpansion::v_SetUpPhysNormals(const int edge)
999  {
1000  boost::ignore_unused(edge);
1001  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
1002  }
1003 
1004  int StdExpansion::v_CalcNumberOfCoefficients(const std::vector<unsigned int> &nummodes, int &modes_offset)
1005  {
1006  boost::ignore_unused(nummodes, modes_offset);
1007  NEKERROR(ErrorUtil::efatal, "This function is not defined for this class");
1008  return 0;
1009  }
1010 
1012  {
1013  boost::ignore_unused(Fx, outarray);
1014  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
1015  }
1016 
1018  {
1019  boost::ignore_unused(Fx, Fy, outarray);
1020  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
1021  }
1022 
1024  const Array<OneD, const NekDouble> &Fy,
1025  const Array<OneD, const NekDouble> &Fz,
1026  Array< OneD, NekDouble> &outarray)
1027  {
1028  boost::ignore_unused(Fx, Fy, Fz, outarray);
1029  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
1030  }
1031 
1033  {
1034  boost::ignore_unused(Fvec, outarray);
1035  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
1036  }
1037 
1038 
1040  {
1041  boost::ignore_unused(mkey);
1042  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
1044  }
1045 
1047  {
1048  boost::ignore_unused(mkey);
1049  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
1050  }
1051 
1053 
1054  {
1055  boost::ignore_unused(face);
1056  NEKERROR(ErrorUtil::efatal, "This function is only valid for three-dimensional LocalRegions");
1057  return eDir1FwdDir1_Dir2FwdDir2;
1058  }
1059 
1061  {
1062  boost::ignore_unused(edge);
1063  NEKERROR(ErrorUtil::efatal, "This function is only valid for two-dimensional LocalRegions");
1064  return eForwards;
1065  }
1066 
1069  Array<OneD, NekDouble> &outarray)
1070  {
1071  boost::ignore_unused(dir, inarray, outarray);
1072  NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1073  }
1074 
1076  Array<OneD, NekDouble> &coeffs,
1078  {
1079  boost::ignore_unused(coeffs, dir);
1080  NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1081  }
1082 
1083 
1085  const Array<OneD, const NekDouble> &physvals)
1086 
1087  {
1088  boost::ignore_unused(Lcoord, physvals);
1089  NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1090  return 0;
1091  }
1092 
1094  {
1095  boost::ignore_unused(xi, eta);
1096  NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1097  }
1098 
1100  {
1101  ASSERTL0(false, "This function is needs defining for this shape");
1102  return 0;
1103  }
1104 
1106  {
1107  ASSERTL0(false, "This function is needs defining for this shape");
1108  return 0;
1109  }
1110 
1112  {
1113  ASSERTL0(false, "This function is needs defining for this shape");
1114  return 0;
1115  }
1116 
1118  {
1119  ASSERTL0(false, "This function is needs defining for this shape");
1120  return 0;
1121  }
1122 
1123  int StdExpansion::v_GetEdgeNcoeffs(const int i) const
1124  {
1125  boost::ignore_unused(i);
1126  ASSERTL0(false, "This function is not valid or not defined");
1127  return 0;
1128  }
1129 
1131  {
1132  ASSERTL0(false, "This function is not valid or not defined");
1133  return 0;
1134  }
1135 
1136  int StdExpansion::v_GetEdgeNumPoints(const int i) const
1137  {
1138  boost::ignore_unused(i);
1139  ASSERTL0(false, "This function is not valid or not defined");
1140  return 0;
1141  }
1142 
1144  {
1145  boost::ignore_unused(edge);
1146  ASSERTL0(false, "This function is not valid or not defined");
1147  return 0;
1148  }
1149 
1151  {
1152  boost::ignore_unused(i);
1153  ASSERTL0(false, "This function is not valid or not defined");
1155  }
1156 
1157  const LibUtilities::BasisKey StdExpansion::v_DetFaceBasisKey(const int i, const int k) const
1158  {
1159  boost::ignore_unused(i, k);
1160  ASSERTL0(false, "This function is not valid or not defined");
1162  }
1163 
1164  int StdExpansion::v_GetFaceNumPoints(const int i) const
1165  {
1166  boost::ignore_unused(i);
1167  ASSERTL0(false, "This function is not valid or not defined");
1168  return 0;
1169  }
1170 
1171  int StdExpansion::v_GetFaceNcoeffs(const int i) const
1172  {
1173  boost::ignore_unused(i);
1174  ASSERTL0(false, "This function is not valid or not defined");
1175  return 0;
1176  }
1177 
1178  int StdExpansion::v_GetFaceIntNcoeffs(const int i) const
1179  {
1180  boost::ignore_unused(i);
1181  ASSERTL0(false, "This function is not valid or not defined");
1182  return 0;
1183  }
1184 
1186  {
1187  ASSERTL0(false, "This function is not valid or not defined");
1188  return 0;
1189  }
1190 
1191  int StdExpansion::v_GetTraceNcoeffs(const int i) const
1192  {
1193  boost::ignore_unused(i);
1194  ASSERTL0(false, "This function is not valid or not defined");
1195  return 0;
1196  }
1197 
1198 
1200  {
1201  boost::ignore_unused(i, j);
1202  ASSERTL0(false, "This function is not valid or not defined");
1204  }
1205 
1207  {
1208  boost::ignore_unused(i);
1209  ASSERTL0(false, "This function is not valid or not defined");
1210 
1212  }
1213 
1215  {
1216  ASSERTL0(false, "This function is not valid or not defined");
1217 
1219  }
1220 
1222  {
1223  ASSERTL0(false, "This expansion does not have a shape type defined");
1225  }
1226 
1227  std::shared_ptr<StdExpansion>
1229  {
1230  ASSERTL0(false,"This method is not defined for this expansion");
1231  StdExpansionSharedPtr returnval;
1232  return returnval;
1233  }
1234 
1235  std::shared_ptr<StdExpansion>
1237  {
1238  ASSERTL0(false,"This method is not defined for this expansion");
1239  StdExpansionSharedPtr returnval;
1240  return returnval;
1241  }
1242 
1244  {
1245  ASSERTL0(false, "This function is not valid or not defined");
1246  return 0;
1247  }
1248 
1250  {
1251  ASSERTL0(false,"This function has not been defined for this expansion");
1252  return false;
1253  }
1254 
1255 
1257  {
1258  return false;
1259  }
1260 
1262  const Array<OneD, const NekDouble>& inarray,
1263  Array<OneD, NekDouble> &outarray)
1264  {
1265  boost::ignore_unused(dir, inarray, outarray);
1266  NEKERROR(ErrorUtil::efatal, "This method has not been defined");
1267  }
1268 
1269 
1270  /**
1271  *
1272  */
1274  const Array<OneD, const NekDouble>& direction,
1275  const Array<OneD, const NekDouble>& inarray,
1276  Array<OneD, NekDouble> &outarray)
1277  {
1278  boost::ignore_unused(direction, inarray, outarray);
1279  NEKERROR(ErrorUtil::efatal, "This method has not been defined");
1280  }
1281 
1282 
1283  /**
1284  *
1285  */
1287  Array<OneD, NekDouble> &outarray)
1288  {
1289  boost::ignore_unused(inarray, outarray);
1290  NEKERROR(ErrorUtil::efatal, "This method has not been defined");
1291  }
1292 
1293 
1294  /**
1295  * @brief Integrates the specified function over the domain.
1296  * @see StdRegions#StdExpansion#Integral.
1297  */
1299  {
1300  boost::ignore_unused(inarray);
1301  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1302  "local expansions");
1303  return 0;
1304  }
1305 
1306 
1307  void StdExpansion::v_AddRobinMassMatrix(const int edgeid, const Array<OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
1308  {
1309  boost::ignore_unused(edgeid, primCoeffs, inoutmat);
1310  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1311  "specific element types");
1312  }
1313 
1315  {
1316  boost::ignore_unused(edgeid, primCoeffs, coeffs);
1317  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1318  "specific element types");
1319  }
1320 
1321  /**
1322  * @brief Calculate the derivative of the physical points
1323  * @see StdRegions#StdExpansion#PhysDeriv
1324  */
1326  Array<OneD, NekDouble> &out_d1,
1327  Array<OneD, NekDouble> &out_d2,
1328  Array<OneD, NekDouble> &out_d3)
1329  {
1330  boost::ignore_unused(inarray, out_d1, out_d2, out_d3);
1331  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1332  "local expansions");
1333  }
1334 
1336  Array<OneD, NekDouble> &out_ds)
1337  {
1338  boost::ignore_unused(inarray, out_ds);
1339  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1340  "local expansions");
1341  }
1343  Array<OneD, NekDouble>& out_dn)
1344  {
1345  boost::ignore_unused(inarray, out_dn);
1346  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1347  "local expansions");
1348  }
1349 
1350  /**
1351  * @brief Calculate the derivative of the physical points in a
1352  * given direction
1353  * @see StdRegions#StdExpansion#PhysDeriv
1354  */
1355  void StdExpansion::v_PhysDeriv(const int dir,
1356  const Array<OneD, const NekDouble>& inarray,
1357  Array<OneD, NekDouble> &out_d0)
1358 
1359  {
1360  boost::ignore_unused(dir, inarray, out_d0);
1361  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1362  "specific element types");
1363  }
1364 
1365  /**
1366  * @brief Physical derivative along a direction vector.
1367  * @see StdRegions#StdExpansion#PhysDirectionalDeriv
1368  */
1370  const Array<OneD, const NekDouble>& direction,
1371  Array<OneD, NekDouble> &outarray)
1372  {
1373  boost::ignore_unused(inarray, direction, outarray);
1374  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1375  "specific element types");
1376  }
1377 
1379  Array<OneD, NekDouble> &out_d1,
1380  Array<OneD, NekDouble> &out_d2,
1381  Array<OneD, NekDouble> &out_d3)
1382  {
1383  boost::ignore_unused(inarray, out_d1, out_d2, out_d3);
1384  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1385  }
1386 
1387  void StdExpansion::v_StdPhysDeriv (const int dir,
1388  const Array<OneD, const NekDouble>& inarray,
1389  Array<OneD, NekDouble> &outarray)
1390  {
1391  boost::ignore_unused(dir, inarray, outarray);
1392  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1393  }
1394 
1395 
1397  {
1398  boost::ignore_unused(coords, physvals);
1399  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1400  return 0;
1401  }
1402 
1403 
1405  {
1406  boost::ignore_unused(I, physvals);
1407  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1408  return 0;
1409  }
1410 
1411 
1412  void StdExpansion::v_FillMode(const int mode, Array<OneD, NekDouble> &outarray)
1413  {
1414  boost::ignore_unused(mode, outarray);
1415  NEKERROR(ErrorUtil::efatal, "This function has not "
1416  "been defined for this shape");
1417  }
1418 
1420  {
1421  boost::ignore_unused(mkey);
1422  NEKERROR(ErrorUtil::efatal, "This function has not "
1423  "been defined for this element");
1424  DNekMatSharedPtr returnval;
1425  return returnval;
1426  }
1427 
1429  {
1430  boost::ignore_unused(mkey);
1431  NEKERROR(ErrorUtil::efatal, "This function has not "
1432  "been defined for this element");
1433  DNekMatSharedPtr returnval;
1434  return returnval;
1435  }
1436 
1438  Array<OneD, NekDouble> &coords_1,
1439  Array<OneD, NekDouble> &coords_2)
1440  {
1441  boost::ignore_unused(coords_0, coords_1, coords_2);
1442  NEKERROR(ErrorUtil::efatal, "Write coordinate definition method");
1443  }
1444 
1446  Array<OneD, NekDouble> &coord)
1447  {
1448  boost::ignore_unused(Lcoord, coord);
1449  NEKERROR(ErrorUtil::efatal, "Write coordinate definition method");
1450  }
1451 
1453  {
1454  NEKERROR(ErrorUtil::efatal, "Write method");
1455  return -1;
1456  }
1457 
1459  {
1460  boost::ignore_unused(outarray);
1461  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1462  }
1463 
1465  {
1466  boost::ignore_unused(outarray);
1467  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1468  }
1469 
1470  int StdExpansion::v_GetVertexMap(const int localVertexId,
1471  bool useCoeffPacking)
1472  {
1473  boost::ignore_unused(localVertexId, useCoeffPacking);
1474  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1475  return 0;
1476  }
1477 
1478  void StdExpansion::v_GetEdgeInteriorMap(const int eid, const Orientation edgeOrient,
1479  Array<OneD, unsigned int> &maparray,
1480  Array<OneD, int> &signarray)
1481  {
1482  boost::ignore_unused(eid, edgeOrient, maparray, signarray);
1483  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1484  }
1485 
1487  const int fid,
1488  const Orientation faceOrient,
1489  int &numModes0,
1490  int &numModes1)
1491  {
1492  boost::ignore_unused(fid, faceOrient, numModes0, numModes1);
1493  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1494  }
1495 
1496  void StdExpansion::v_GetFaceInteriorMap(const int fid, const Orientation faceOrient,
1497  Array<OneD, unsigned int> &maparray,
1498  Array<OneD, int> &signarray)
1499  {
1500  boost::ignore_unused(fid, faceOrient, maparray, signarray);
1501  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1502  }
1503 
1505  const int eid,
1506  const Orientation edgeOrient,
1507  Array<OneD, unsigned int>& maparray,
1508  Array<OneD, int>& signarray,
1509  int P)
1510  {
1511  boost::ignore_unused(eid, edgeOrient, maparray, signarray, P);
1512  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1513  }
1514 
1515  void StdExpansion::v_GetFaceToElementMap(const int fid, const Orientation faceOrient,
1516  Array<OneD, unsigned int> &maparray,
1517  Array<OneD, int> &signarray,
1518  int nummodesA, int nummodesB)
1519  {
1520  boost::ignore_unused(fid, faceOrient, maparray, signarray,
1521  nummodesA, nummodesB);
1522  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1523  }
1524 
1526  {
1527  boost::ignore_unused(edge, inarray, outarray);
1528  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape or library" );
1529  }
1530 
1531  void StdExpansion::v_GetEdgePhysVals(const int edge, const std::shared_ptr<StdExpansion> &EdgeExp, const Array<OneD, const NekDouble> &inarray, Array<OneD,NekDouble> &outarray)
1532  {
1533  boost::ignore_unused(edge, EdgeExp, inarray, outarray);
1534  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape or library" );
1535  }
1536 
1537  void StdExpansion::v_GetTracePhysVals(const int edge, const std::shared_ptr<StdExpansion> &EdgeExp, const Array<OneD, const NekDouble> &inarray, Array<OneD,NekDouble> &outarray, StdRegions::Orientation orient)
1538  {
1539  boost::ignore_unused(edge, EdgeExp, inarray, outarray, orient);
1540  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape or library" );
1541  }
1542 
1543  void StdExpansion::v_GetVertexPhysVals(const int vertex, const Array<OneD, const NekDouble> &inarray, NekDouble &outarray)
1544  {
1545  boost::ignore_unused(vertex, inarray, outarray);
1546  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape or library" );
1547  }
1548 
1550  {
1551  boost::ignore_unused(edge, inarray, outarray);
1552  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape or library" );
1553  }
1554 
1556  const int edge,
1557  Array<OneD, NekDouble> &outarray)
1558  {
1559  boost::ignore_unused(edge, outarray);
1561  "Method does not exist for this shape or library");
1562  }
1563 
1564  void StdExpansion::v_GetFacePhysVals( const int face,
1565  const std::shared_ptr<StdExpansion> &FaceExp,
1566  const Array<OneD, const NekDouble> &inarray,
1567  Array<OneD, NekDouble> &outarray,
1568  StdRegions::Orientation orient)
1569  {
1570  boost::ignore_unused(face, FaceExp, inarray, outarray, orient);
1571  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape or library" );
1572  }
1573 
1575  const int edge,
1576  Array<OneD, int> &outarray)
1577  {
1578  boost::ignore_unused(edge, outarray);
1580  "Method does not exist for this shape or library" );
1581  }
1582 
1583  void StdExpansion::v_GetFacePhysMap(const int face,
1584  Array<OneD, int> &outarray)
1585  {
1586  boost::ignore_unused(face, outarray);
1587  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape or library" );
1588  }
1589 
1591  const Array<OneD, const NekDouble> &inarray,
1592  Array<OneD, NekDouble> &outarray)
1593  {
1594  boost::ignore_unused(inarray, outarray);
1595  v_MultiplyByStdQuadratureMetric(inarray,outarray);
1596  }
1597 
1599  const Array<OneD, const NekDouble> &inarray,
1600  Array<OneD, NekDouble> &outarray)
1601  {
1602  boost::ignore_unused(inarray, outarray);
1603  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape or library");
1604  }
1605 
1607  Array<OneD, NekDouble> &outarray)
1608  {
1609  boost::ignore_unused(inarray, outarray);
1610  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1611  }
1612 
1614  Array<OneD, NekDouble> &outarray,
1615  bool multiplybyweights)
1616  {
1617  boost::ignore_unused(inarray, outarray, multiplybyweights);
1618  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1619  }
1620 
1621  /**
1622  *
1623  */
1625  const Array<OneD, const NekDouble>& direction,
1626  const Array<OneD, const NekDouble>& inarray,
1627  Array<OneD, NekDouble> &outarray)
1628  {
1629  boost::ignore_unused(direction, inarray, outarray);
1631  "Method does not exist for this shape" );
1632  }
1633 
1635  const Array<OneD, const NekDouble>& inarray,
1636  Array<OneD, NekDouble> &outarray)
1637  {
1638  boost::ignore_unused(dir, inarray, outarray);
1639  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1640  }
1641 
1643  Array<OneD,NekDouble> &outarray,
1644  const StdMatrixKey &mkey)
1645  {
1646  // If this function is not reimplemented on shape level, the function
1647  // below will be called
1648  MassMatrixOp_MatFree(inarray,outarray,mkey);
1649  }
1650 
1652  Array<OneD,NekDouble> &outarray,
1653  const StdMatrixKey &mkey)
1654  {
1655  // If this function is not reimplemented on shape level, the function
1656  // below will be called
1657  LaplacianMatrixOp_MatFree(inarray,outarray,mkey);
1658  }
1659 
1660 
1662  const StdMatrixKey &mkey)
1663  {
1664  boost::ignore_unused(array, mkey);
1665  ASSERTL0(false, "This function is not defined in StdExpansion.");
1666  }
1667 
1669  Array<OneD, NekDouble> &array,
1670  const NekDouble alpha,
1671  const NekDouble exponent,
1672  const NekDouble cutoff)
1673  {
1674  boost::ignore_unused(array, alpha, exponent, cutoff);
1675  ASSERTL0(false, "This function is not defined in StdExpansion.");
1676  }
1677 
1679  const Array<OneD, const NekDouble> &inarray,
1680  Array<OneD, NekDouble> &outarray)
1681  {
1682  boost::ignore_unused(numMin, inarray, outarray);
1683  ASSERTL0(false, "This function is not defined in StdExpansion.");
1684  }
1685 
1686  void StdExpansion::v_LaplacianMatrixOp(const int k1, const int k2,
1687  const Array<OneD, const NekDouble> &inarray,
1688  Array<OneD,NekDouble> &outarray,
1689  const StdMatrixKey &mkey)
1690  {
1691  // If this function is not reimplemented on shape level, the function
1692  // below will be called
1693  LaplacianMatrixOp_MatFree(k1,k2,inarray,outarray,mkey);
1694  }
1695 
1697  const Array<OneD, const NekDouble> &inarray,
1698  Array<OneD,NekDouble> &outarray,
1699  const StdMatrixKey &mkey)
1700  {
1701  // If this function is not reimplemented on shape level, the function
1702  // below will be called
1703  WeakDerivMatrixOp_MatFree(i,inarray,outarray,mkey);
1704 
1705  }
1706 
1708  Array<OneD,NekDouble> &outarray,
1709  const StdMatrixKey &mkey)
1710  {
1711  // If this function is not reimplemented on shape level, the function
1712  // below will be called
1713  WeakDirectionalDerivMatrixOp_MatFree(inarray,outarray,mkey);
1714 
1715  }
1716 
1718  Array<OneD,NekDouble> &outarray,
1719  const StdMatrixKey &mkey)
1720  {
1721  // If this function is not reimplemented on shape level, the function
1722  // below will be called
1723  MassLevelCurvatureMatrixOp_MatFree(inarray,outarray,mkey);
1724  }
1725 
1727  const NekDouble> &inarray,
1728  Array<OneD,NekDouble> &outarray,
1729  const StdMatrixKey &mkey, bool addDiffusionTerm)
1730  {
1731  // If this function is not reimplemented on shape level, the function
1732  // below will be called
1733  LinearAdvectionDiffusionReactionMatrixOp_MatFree(inarray,outarray,mkey,addDiffusionTerm);
1734 
1735  }
1736 
1738  Array<OneD,NekDouble> &outarray,
1739  const StdMatrixKey &mkey)
1740  {
1741  // If this function is not reimplemented on shape level, the function
1742  // below will be called
1743  HelmholtzMatrixOp_MatFree(inarray,outarray,mkey);
1744  }
1745 
1747  Array<OneD,NekDouble> &outarray,
1748  const StdMatrixKey &mkey)
1749  {
1750  // If this function is not reimplemented on shape level, the function
1751  // below will be called
1752  LaplacianMatrixOp_MatFree_GenericImpl(inarray,outarray,mkey);
1753  }
1754 
1756  const Array<OneD, const NekDouble> &inarray,
1757  Array<OneD, NekDouble> &outarray,
1759  {
1760  boost::ignore_unused(inarray, outarray, wsp);
1761  ASSERTL0(false, "Not implemented.");
1762  }
1763 
1765  Array<OneD,NekDouble> &outarray,
1766  const StdMatrixKey &mkey)
1767  {
1768  // If this function is not reimplemented on shape level, the function
1769  // below will be called
1770  HelmholtzMatrixOp_MatFree_GenericImpl(inarray,outarray,mkey);
1771  }
1772 
1773  const NormalVector & StdExpansion::v_GetEdgeNormal(const int edge) const
1774  {
1775  boost::ignore_unused(edge);
1776  ASSERTL0(false, "Cannot get edge normals for this expansion.");
1777  static NormalVector result;
1778  return result;
1779  }
1780 
1782  {
1783  boost::ignore_unused(edge);
1784  ASSERTL0(false, "Cannot compute edge normal for this expansion.");
1785  }
1786 
1788  {
1789  boost::ignore_unused(edge);
1790  ASSERTL0(false, "Not implemented.");
1791  }
1792 
1794  {
1795  boost::ignore_unused(edge);
1796  ASSERTL0(false, "Not implemented.");
1797  return false;
1798  }
1799 
1801  {
1802  boost::ignore_unused(face);
1803  ASSERTL0(false, "Cannot compute face normal for this expansion.");
1804  }
1805 
1807  {
1808  boost::ignore_unused(face);
1809  ASSERTL0(false, "Not implemented.");
1810  }
1811 
1813  {
1814  boost::ignore_unused(face);
1815  ASSERTL0(false, "Not implemented.");
1816  return false;
1817  }
1818 
1820  {
1821  boost::ignore_unused(vertex);
1822  ASSERTL0(false, "Cannot compute vertex normal for this expansion.");
1823  }
1824 
1825  void StdExpansion::v_NegateVertexNormal(const int vertex)
1826  {
1827  boost::ignore_unused(vertex);
1828  ASSERTL0(false, "Not implemented.");
1829  }
1830 
1832  {
1833  boost::ignore_unused(vertex);
1834  ASSERTL0(false, "Not implemented.");
1835  return false;
1836  }
1837 
1838  const NormalVector & StdExpansion::v_GetFaceNormal(const int face) const
1839  {
1840  boost::ignore_unused(face);
1841  ASSERTL0(false, "Cannot get face normals for this expansion.");
1842  static NormalVector result;
1843  return result;
1844  }
1845 
1846  const NormalVector & StdExpansion::v_GetVertexNormal(const int vertex) const
1847  {
1848  boost::ignore_unused(vertex);
1849  ASSERTL0(false, "Cannot get vertex normals for this expansion.");
1850  static NormalVector result;
1851  return result;
1852  }
1853 
1855  {
1856  boost::ignore_unused(id);
1857  ASSERTL0(false, "Cannot get face normals for this expansion.");
1858  static NormalVector result;
1859  return result;
1860  }
1861 
1864  {
1865  boost::ignore_unused(eid);
1866  ASSERTL0(false, "Not implemented.");
1867  Array<OneD, unsigned int> noinversemap(1);
1868  return noinversemap;
1869  }
1870 
1873  StdRegions::Orientation faceOrient,
1874  int P1,
1875  int P2)
1876  {
1877  boost::ignore_unused(fid, faceOrient, P1, P2);
1878  ASSERTL0(false, "Not implemented.");
1879  Array<OneD, unsigned int> noinversemap(1);
1880  return noinversemap;
1881  }
1882 
1887  {
1888  boost::ignore_unused(vmap, emap, fmap);
1889  ASSERTL0(false, "Not implemented.");
1890  }
1891 
1894  const DNekScalMatSharedPtr & m_transformationmatrix)
1895  {
1896  boost::ignore_unused(m_transformationmatrix);
1897  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
1898  return NullDNekMatSharedPtr;
1899  }
1900 
1902  const Array<OneD, const NekDouble> &inarray,
1903  Array<OneD, NekDouble> &outarray,
1904  int npset)
1905  {
1907  DNekMatSharedPtr intmat;
1908 
1909  int nqtot = GetTotPoints();
1910  int np = 0;
1911  if(npset == -1) // use values from basis num points()
1912  {
1913  int nqbase;
1914  for(int i = 0; i < m_base.num_elements(); ++i)
1915  {
1916  nqbase = m_base[i]->GetNumPoints();
1917  np = std::max(np,nqbase);
1918  }
1919 
1920  StdMatrixKey Ikey(ePhysInterpToEquiSpaced, shape, *this);
1921  intmat = GetStdMatrix(Ikey);
1922  }
1923  else
1924  {
1925  np = npset;
1926 
1927  ConstFactorMap cmap;
1928  cmap[eFactorConst] = np;
1929  StdMatrixKey Ikey(ePhysInterpToEquiSpaced, shape, *this, cmap);
1930  intmat = GetStdMatrix(Ikey);
1931 
1932  }
1933 
1934  NekVector<NekDouble> in (nqtot,inarray,eWrapper);
1936  out = (*intmat) * in;
1937  }
1938 
1940  Array<OneD, int> &conn,
1941  bool standard)
1942  {
1943  boost::ignore_unused(conn, standard);
1944  ASSERTL0(false, "Not implemented.");
1945  }
1946 
1948  const Array<OneD, const NekDouble> &inarray,
1949  Array<OneD, NekDouble> &outarray)
1950  {
1952 
1953  // inarray has to be consistent with NumModes definition
1954  // There is also a check in GetStdMatrix to see if all
1955  // modes are of the same size
1956  ConstFactorMap cmap;
1957 
1958  cmap[eFactorConst] = m_base[0]->GetNumModes();
1959  StdMatrixKey Ikey(eEquiSpacedToCoeffs, shape, *this,cmap);
1960  DNekMatSharedPtr intmat = GetStdMatrix(Ikey);
1961 
1962  NekVector<NekDouble> in (m_ncoeffs, inarray, eWrapper);
1963  NekVector<NekDouble> out(m_ncoeffs, outarray,eWrapper);
1964  out = (*intmat) * in;
1965  }
1966 
1967  }//end namespace
1968 }//end namespace
virtual bool v_EdgeNormalNegated(const int edge)
void HelmholtzMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_GetFaceNumModes(const int fid, const Orientation faceOrient, int &numModes0, int &numModes1)
virtual const LibUtilities::BasisKey v_DetFaceBasisKey(const int i, const int k) const
static const PointsKey NullPointsKey(0, eNoPointsType)
virtual int v_GetVertexMap(int localVertexId, bool useCoeffPacking=false)
virtual void v_ComputeFaceNormal(const int face)
virtual LibUtilities::BasisType v_GetEdgeBasisType(const int i) const
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
virtual void v_GetInteriorMap(Array< OneD, unsigned int > &outarray)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:209
virtual int v_GetTotalFaceIntNcoeffs() const
void LaplacianMatrixOp_MatFree_GenericImpl(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)
virtual void v_NegateFaceNormal(const int face)
DNekMatSharedPtr CreateStdMatrix(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:597
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:974
static Array< OneD, NekDouble > NullNekDouble1DArray
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:469
virtual void v_FwdTrans_BndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:16
void IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
virtual void v_GetEdgePhysVals(const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Extract the physical values along edge edge from inarray into outarray following the local edge orien...
virtual void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray)
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:73
virtual void v_GetFaceInteriorMap(const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
virtual int v_NumBndryCoeffs() const
virtual int v_GetFaceNumPoints(const int i) const
void WeakDirectionalDerivMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_PhysDeriv_n(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn)
virtual void v_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
Definition: StdExpansion.h:995
LibUtilities::NekManager< IndexMapKey, IndexMapValues, IndexMapKey::opLess > m_IndexMapManager
virtual int v_GetShapeDimension() const
virtual void v_MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray)
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:488
virtual int v_GetTraceNcoeffs(const int i) const
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:445
Orientation GetIndexOrientation() const
Definition: IndexMapKey.h:91
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual const NormalVector & v_GetFaceNormal(const int face) const
STL namespace.
static const NekDouble kNekSqrtTol
virtual void v_DropLocStaticCondMatrix(const LocalRegions::MatrixKey &mkey)
virtual ~StdExpansion()
Destructor.
virtual void v_GetEdgeToElementMap(const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int P=-1)
virtual int v_GetFaceIntNcoeffs(const int i) const
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:294
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey)
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.
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
LibUtilities::NekManager< StdMatrixKey, DNekBlkMat, StdMatrixKey::opLess > m_stdStaticCondMatrixManager
virtual StdRegions::Orientation v_GetForient(int face)
virtual void v_GetEdgeInteriorMap(const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:81
const Array< OneD, const NekDouble > & GetVarCoeff(const StdRegions::VarCoeffType &coeff) const
Definition: StdMatrixKey.h:145
virtual void v_NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
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.
virtual int v_GetFaceNcoeffs(const int i) const
void Vabs(int n, const T *x, const int incx, T *y, const int incy)
vabs: y = |x|
Definition: Vmath.cpp:427
bool HasVarCoeff(const StdRegions::VarCoeffType &coeff) const
Definition: StdMatrixKey.h:166
virtual const NormalVector & v_GetSurfaceNormal(const int id) const
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:71
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:714
virtual int v_GetEdgeNcoeffs(const int i) const
void 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_GetTracePhysVals(const int edge, const std::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=eNoOrientation)
int GetNumberOfCoefficients(ShapeType shape, std::vector< unsigned int > &modes, int offset)
Definition: ShapeType.hpp:313
virtual void v_GetEdgePhysMap(const int edge, Array< OneD, int > &outarray)
static DNekMatSharedPtr NullDNekMatSharedPtr
Definition: NekTypeDefs.hpp:78
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 bool v_VertexNormalNegated(const int vertex)
virtual void v_AddRobinMassMatrix(const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
void WeakDerivMatrixOp(const int i, 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)
BasisManagerT & BasisManager(void)
virtual void v_GetInverseBoundaryMaps(Array< OneD, unsigned int > &vmap, Array< OneD, Array< OneD, unsigned int > > &emap, Array< OneD, Array< OneD, unsigned int > > &fmap)
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_MassLevelCurvatureMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LinearAdvectionDiffusionReactionMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual int v_GetElmtId()
Get the element id of this expansion when used in a list by returning value of m_elmt_id.
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
virtual void v_SetUpPhysNormals(const int edge)
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray)
Integrates the specified function over the domain.
virtual void v_GetEdgeInterpVals(const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
The base class for all shapes.
Definition: StdExpansion.h:68
virtual void SetUpPhysNormals(const int edge)
virtual int v_GetEdgeNumPoints(const int i) const
virtual void v_NegateEdgeNormal(const int edge)
virtual void v_AddRobinEdgeContribution(const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, Array< OneD, NekDouble > &coeffs)
void GeneralMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_NegateVertexNormal(const int vertex)
virtual void v_GetCoord(const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey)
virtual LibUtilities::ShapeType v_DetShapeType() const
virtual std::shared_ptr< StdExpansion > v_GetLinStdExp(void) const
IndexMapValuesSharedPtr CreateIndexMap(const IndexMapKey &ikey)
Create an IndexMap which contains mapping information linking any specific element shape with either ...
IndexMapType GetIndexMapType() const
Definition: IndexMapKey.h:86
virtual void v_PhysDeriv_s(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds)
virtual void v_GetSimplexEquiSpacedConnectivity(Array< OneD, int > &conn, bool standard=true)
virtual const NormalVector & v_GetEdgeNormal(const int edge) const
virtual void v_WeakDirectionalDerivMatrixOp(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 NekDouble v_PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
virtual void v_GetVertexPhysVals(const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray)
void WeakDerivMatrixOp_MatFree(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:817
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
virtual void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
Defines a specification for a set of points.
Definition: Points.h:59
virtual const LibUtilities::BasisKey v_DetEdgeBasisKey(const int i) const
double NekDouble
virtual bool v_IsBoundaryInteriorExpansion()
static void Dgemv(const char &trans, const int &m, const int &n, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = A x where A[m x n].
Definition: Blas.hpp:168
void BwdTrans_MatOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual LibUtilities::PointsKey v_GetFacePointsKey(const int i, const int j) const
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:828
virtual void v_GetEdgeQFactors(const int edge, Array< OneD, NekDouble > &outarray)
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_StdPhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2, Array< OneD, NekDouble > &out_d3)
virtual void v_GetFaceToElementMap(const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int nummodesA=-1, int nummodesB=-1)
NekDouble StdPhysEvaluate(const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
NekDouble GetConstFactor(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:121
virtual void v_LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
void HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
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.
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:130
virtual int v_NumDGBndryCoeffs() const
void MassLevelCurvatureMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual Array< OneD, unsigned int > v_GetEdgeInverseBoundaryMap(int eid)
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix(const LocalRegions::MatrixKey &mkey)
virtual void v_ComputeEdgeNormal(const int edge)
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:346
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...
virtual int v_DetCartesianDirOfEdge(const int edge)
void LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:981
virtual void v_SetCoeffsToOrientation(StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual const LibUtilities::PointsKey v_GetNodalPointsKey() const
virtual void v_ComputeVertexNormal(const int vertex)
virtual void v_IProductWRTDirectionalDerivBase_SumFac(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_GetFacePhysMap(const int face, Array< OneD, int > &outarray)
virtual NekDouble v_StdPhysEvaluate(const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
bool ConstFactorExists(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:130
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:274
virtual StdRegions::Orientation v_GetEorient(int edge)
void MassMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
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.
virtual int v_CalcNumberOfCoefficients(const std::vector< unsigned int > &nummodes, int &modes_offset)
void HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_WeakDerivMatrixOp(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_GetFacePhysVals(const int face, const std::shared_ptr< StdExpansion > &FaceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient)
virtual bool v_FaceNormalNegated(const int face)
virtual int v_GetTotalEdgeIntNcoeffs() const
virtual void v_SetPhysNormals(Array< OneD, const NekDouble > &normal)
DNekBlkMatSharedPtr CreateStdStaticCondMatrix(const StdMatrixKey &mkey)
Create the static condensation of a matrix when using a boundary interior decomposition.
static DNekScalBlkMatSharedPtr NullDNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:80
virtual void v_ExponentialFilter(Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.
virtual const Array< OneD, const NekDouble > & v_GetPhysNormals(void)
LibUtilities::NekManager< StdMatrixKey, DNekMat, StdMatrixKey::opLess > m_stdMatrixManager
virtual const NormalVector & v_GetVertexNormal(const int vertex) const
StdExpansion()
Default Constructor.
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140
std::shared_ptr< IndexMapValues > IndexMapValuesSharedPtr
Definition: IndexMapKey.h:124
void WeakDirectionalDerivMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:376
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...
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
Array< OneD, LibUtilities::BasisSharedPtr > m_base
virtual Array< OneD, unsigned int > v_GetFaceInverseBoundaryMap(int fid, StdRegions::Orientation faceOrient=eNoOrientation, int P1=-1, int P2=-1)
void GeneralMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
LibUtilities::PointsType GetNodalPointsType() const
Definition: StdMatrixKey.h:91
std::shared_ptr< StdMatrixKey > StdMatrixKeySharedPtr
Definition: StdMatrixKey.h:189
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:110
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
virtual void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
LibUtilities::ShapeType GetShapeType() const
Definition: StdMatrixKey.h:86
virtual void v_LinearAdvectionDiffusionReactionMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_0, Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2)
virtual void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual std::shared_ptr< StdExpansion > v_GetStdExp(void) const
void GetBoundaryMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:812
virtual DNekMatSharedPtr v_BuildInverseTransformationMatrix(const DNekScalMatSharedPtr &m_transformationmatrix)
Describes the specification for a Basis.
Definition: Basis.h:49
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)=0
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:302
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:186
void BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:265
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:295
void LinearAdvectionDiffusionReactionMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)