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 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Definition of methods in class StdExpansion which is
33 // the base class to all expansion shapes
34 //
35 ///////////////////////////////////////////////////////////////////////////////
36 
38 #include <LibUtilities/Foundations/ManagerAccess.h> // for BasisManager, etc
39 
40 namespace Nektar
41 {
42  namespace StdRegions
43  {
45  m_elmt_id(0),
46  m_ncoeffs(0)
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),
55  m_elmt_id(0),
56  m_ncoeffs(numcoeffs),
57  m_stdMatrixManager(
58  boost::bind(&StdExpansion::CreateStdMatrix, this, _1),
59  std::string("StdExpansionStdMatrix")),
60  m_stdStaticCondMatrixManager(
61  boost::bind(&StdExpansion::CreateStdStaticCondMatrix, this, _1),
62  std::string("StdExpansionStdStaticCondMatrix")),
63  m_IndexMapManager(
64  boost::bind(&StdExpansion::CreateIndexMap,this, _1),
65  std::string("StdExpansionIndexMap"))
66  {
67  switch(m_base.num_elements())
68  {
69  case 3:
71  "NULL Basis attempting to be used.");
73 
74  case 2:
76  "NULL Basis attempting to be used.");
77 
79  case 1:
81  "NULL Basis attempting to be used.");
83  break;
84  default:
85  break;
86 // ASSERTL0(false, "numbases incorrectly specified");
87  };
88 
89  } //end constructor
90 
91 
93  m_base(T.m_base),
94  m_elmt_id(T.m_elmt_id),
95  m_ncoeffs(T.m_ncoeffs),
96  m_stdMatrixManager(T.m_stdMatrixManager),
97  m_stdStaticCondMatrixManager(T.m_stdStaticCondMatrixManager),
98  m_IndexMapManager(T.m_IndexMapManager)
99  {
100  }
101 
103  {
104  }
105 
107  const Array<OneD, const NekDouble>& sol)
108  {
109  NekDouble val;
110  int ntot = GetTotPoints();
111  Array<OneD, NekDouble> wsp(ntot);
112 
113  if(sol == NullNekDouble1DArray)
114  {
115  Vmath::Vabs(ntot, phys, 1, wsp, 1);
116  }
117  else
118  {
119  Vmath::Vsub(ntot, sol, 1, phys, 1, wsp, 1);
120  Vmath::Vabs(ntot, wsp, 1, wsp, 1);
121  }
122 
123  val = Vmath::Vamax(ntot, wsp, 1);
124 
125  return val;
126  }
127 
129  const Array<OneD, const NekDouble>& sol)
130  {
131  NekDouble val;
132  int ntot = GetTotPoints();
133  Array<OneD, NekDouble> wsp(ntot);
134 
135  if (sol.num_elements() == 0)
136  {
137  Vmath::Vmul(ntot, phys, 1, phys, 1, wsp, 1);
138  }
139  else
140  {
141  Vmath::Vsub(ntot, sol, 1, phys, 1, wsp, 1);
142  Vmath::Vmul(ntot, wsp, 1, wsp, 1, wsp, 1);
143  }
144 
145  val = v_Integral(wsp);
146 
147  // if val too small, sqrt returns nan.
149  {
150  return 0.0;
151  }
152  else
153  {
154  return sqrt(val);
155  }
156  }
157 
159  const Array<OneD, const NekDouble>& sol)
160  {
161  int i;
162  NekDouble val;
163  int ntot = GetTotPoints();
164  int coordim = v_GetCoordim();
165  Array<OneD, NekDouble> wsp(3*ntot);
166  Array<OneD, NekDouble> wsp_deriv = wsp + ntot;
167  Array<OneD, NekDouble> sum = wsp_deriv + ntot;
168 
169  if(sol == NullNekDouble1DArray)
170  {
171  Vmath::Vcopy(ntot,phys, 1, wsp, 1);
172  Vmath::Vmul(ntot, phys, 1, phys, 1, sum, 1);
173  }
174  else
175  {
176  Vmath::Vsub(ntot, sol, 1, phys, 1, wsp, 1);
177  Vmath::Vmul(ntot, wsp, 1, wsp, 1, sum, 1);
178  }
179 
180 
181  for(i = 0; i < coordim; ++i)
182  {
183  v_PhysDeriv(i,wsp,wsp_deriv);
184  Vmath::Vvtvp(ntot,wsp_deriv,1,wsp_deriv,1,sum,1,sum,1);
185  }
186 
187  val = sqrt(v_Integral(sum));
188 
189  return val;
190  }
191 
192 
194  {
195  DNekBlkMatSharedPtr returnval;
196 
197  DNekMatSharedPtr mat = GetStdMatrix(mkey);
198  int nbdry = NumBndryCoeffs(); // also checks to see if this is a boundary interior decomposed expansion
199  int nint = m_ncoeffs - nbdry;
204 
205  int i,j;
206 
207  Array<OneD,unsigned int> bmap(nbdry);
208  Array<OneD,unsigned int> imap(nint);
209  GetBoundaryMap(bmap);
210  GetInteriorMap(imap);
211 
212  for(i = 0; i < nbdry; ++i)
213  {
214  for(j = 0; j < nbdry; ++j)
215  {
216  (*A)(i,j) = (*mat)(bmap[i],bmap[j]);
217  }
218 
219  for(j = 0; j < nint; ++j)
220  {
221  (*B)(i,j) = (*mat)(bmap[i],imap[j]);
222  }
223  }
224 
225  for(i = 0; i < nint; ++i)
226  {
227  for(j = 0; j < nbdry; ++j)
228  {
229  (*C)(i,j) = (*mat)(imap[i],bmap[j]);
230  }
231 
232  for(j = 0; j < nint; ++j)
233  {
234  (*D)(i,j) = (*mat)(imap[i],imap[j]);
235  }
236  }
237 
238  // Calculate static condensed system
239  if(nint)
240  {
241  D->Invert();
242  (*B) = (*B)*(*D);
243  (*A) = (*A) - (*B)*(*C);
244  }
245 
246  // set up block matrix system
247  Array<OneD, unsigned int> exp_size(2);
248  exp_size[0] = nbdry;
249  exp_size[1] = nint;
250  returnval = MemoryManager<DNekBlkMat>::AllocateSharedPtr(exp_size,exp_size);
251 
252  returnval->SetBlock(0,0,A);
253  returnval->SetBlock(0,1,B);
254  returnval->SetBlock(1,0,C);
255  returnval->SetBlock(1,1,D);
256 
257  return returnval;
258  }
259 
261  {
262  IndexMapValuesSharedPtr returnval;
263 
264  IndexMapType itype = ikey.GetIndexMapType();
265 
266  int entity = ikey.GetIndexEntity();
267 
268  Orientation orient = ikey.GetIndexOrientation();
269 
272 
273  switch(itype)
274  {
275  case eEdgeToElement:
276  {
277  v_GetEdgeToElementMap(entity,orient,map,sign);
278  }
279  break;
280  case eFaceToElement:
281  {
282  v_GetFaceToElementMap(entity,orient,map,sign);
283  }
284  break;
285  case eEdgeInterior:
286  {
287  v_GetEdgeInteriorMap(entity,orient,map,sign);
288  }
289  break;
290  case eFaceInterior:
291  {
292  v_GetFaceInteriorMap(entity,orient,map,sign);
293  }
294  break;
295  case eBoundary:
296  {
297  ASSERTL0(false,"Boundary Index Map not implemented yet.");
298  }
299  break;
300  case eVertex:
301  {
302  ASSERTL0(false,"Vertex Index Map not implemented yet.");
303  }
304  break;
305  default:
306  {
307  ASSERTL0(false,"The Index Map you are requiring is not between the possible options.");
308  }
309  }
310 
311  returnval = MemoryManager<IndexMapValues>::AllocateSharedPtr(map.num_elements());
312 
313  for(int i = 0; i < map.num_elements(); i++)
314  {
315  (*returnval)[i].index = map[i];
316  (*returnval)[i].sign = sign[i];
317  }
318 
319  return returnval;
320  }
321 
323  {
324  int i;
325  DNekMatSharedPtr returnval;
326 
327  switch(mkey.GetMatrixType())
328  {
329  case eInvMass:
330  {
332  DNekMatSharedPtr mmat = GetStdMatrix(masskey);
333 
334  returnval = MemoryManager<DNekMat>::AllocateSharedPtr(*mmat); //Populate standard mass matrix.
335  returnval->Invert();
336  }
337  break;
338  case eInvNBasisTrans:
339  {
341  DNekMatSharedPtr tmpmat = GetStdMatrix(tmpkey);
342  returnval = MemoryManager<DNekMat>::AllocateSharedPtr(*tmpmat); //Populate matrix.
343  returnval->Invert();
344  }
345  break;
346  case eBwdTrans:
347  {
348  int nq = GetTotPoints();
350  Array<OneD, NekDouble> tmpout(nq);
351 
353 
354  for(int i=0; i<m_ncoeffs; ++i)
355  {
356  Vmath::Zero(m_ncoeffs, tmpin, 1);
357  tmpin[i] = 1.0;
358 
359  BwdTrans_SumFac(tmpin,tmpout);
360 
361  Vmath::Vcopy(nq,tmpout.get(),1,
362  returnval->GetRawPtr()+i*nq,1);
363  }
364  }
365  break;
366  case eIProductWRTBase:
367  {
368  int nq = GetTotPoints();
369  Array<OneD, NekDouble> tmpin(nq);
371 
373 
374  for(i=0; i < nq; ++i)
375  {
376  Vmath::Zero(nq, tmpin, 1);
377  tmpin[i] = 1.0;
378 
379  IProductWRTBase_SumFac(tmpin,tmpout);
380 
381  Vmath::Vcopy(m_ncoeffs,tmpout.get(),1,
382  returnval->GetRawPtr()+i*m_ncoeffs,1);
383  }
384  }
385  break;
387  {
388  int nq = GetTotPoints();
389  Array<OneD, NekDouble> tmpin(nq);
391 
393 
394  for(i=0; i < nq; ++i)
395  {
396  Vmath::Zero(nq, tmpin, 1);
397  tmpin[i] = 1.0;
398 
399  IProductWRTDerivBase_SumFac(0,tmpin,tmpout);
400 
401  Vmath::Vcopy(m_ncoeffs,tmpout.get(),1,
402  returnval->GetRawPtr()+i*m_ncoeffs,1);
403  }
404  }
405  break;
407  {
408  int nq = GetTotPoints();
409  Array<OneD, NekDouble> tmpin(nq);
411 
413 
414  for(i=0; i < nq; ++i)
415  {
416  Vmath::Zero(nq, tmpin, 1);
417  tmpin[i] = 1.0;
418 
419  IProductWRTDerivBase_SumFac(1,tmpin,tmpout);
420 
421  Vmath::Vcopy(m_ncoeffs,tmpout.get(),1,
422  returnval->GetRawPtr()+i*m_ncoeffs,1);
423  }
424  }
425  break;
427  {
428  int nq = GetTotPoints();
429  Array<OneD, NekDouble> tmpin(nq);
431 
433 
434  for(i=0; i < nq; ++i)
435  {
436  Vmath::Zero(nq, tmpin, 1);
437  tmpin[i] = 1.0;
438 
439  IProductWRTDerivBase_SumFac(2,tmpin,tmpout);
440 
441  Vmath::Vcopy(m_ncoeffs,tmpout.get(),1,
442  returnval->GetRawPtr()+i*m_ncoeffs,1);
443  }
444  }
445  break;
446  case eMass:
447  case eHelmholtz:
448  case eLaplacian:
449  case eLaplacian00:
450  case eLaplacian01:
451  case eLaplacian02:
452  case eLaplacian11:
453  case eLaplacian12:
454  case eLaplacian22:
455  case eWeakDeriv0:
456  case eWeakDeriv1:
457  case eWeakDeriv2:
459  case eMassLevelCurvature:
462  {
465  DNekMat &Mat = *returnval;
466 
467  for(i=0; i < m_ncoeffs; ++i)
468  {
469  Vmath::Zero(m_ncoeffs, tmp, 1);
470  tmp[i] = 1.0;
471 
472  GeneralMatrixOp_MatFree(tmp,tmp,mkey);
473 
474  Vmath::Vcopy(m_ncoeffs,&tmp[0],1,
475  &(Mat.GetPtr())[0]+i*m_ncoeffs,1);
476  }
477  }
478  break;
479  default:
480  {
481  NEKERROR(ErrorUtil::efatal, "This type of matrix can not be created using a general approach");
482  }
483  break;
484  }
485 
486  return returnval;
487  }
488 
490  Array<OneD,NekDouble> &outarray,
491  const StdMatrixKey &mkey)
492  {
493  switch(mkey.GetMatrixType())
494  {
495  case eMass:
496  MassMatrixOp(inarray,outarray,mkey);
497  break;
498  case eWeakDeriv0:
499  WeakDerivMatrixOp(0,inarray,outarray,mkey);
500  break;
501  case eWeakDeriv1:
502  WeakDerivMatrixOp(1,inarray,outarray,mkey);
503  break;
504  case eWeakDeriv2:
505  WeakDerivMatrixOp(2,inarray,outarray,mkey);
506  break;
508  WeakDirectionalDerivMatrixOp(inarray,outarray,mkey);
509  break;
510  case eMassLevelCurvature:
511  MassLevelCurvatureMatrixOp(inarray,outarray,mkey);
512  break;
514  LinearAdvectionDiffusionReactionMatrixOp(inarray,outarray,mkey,false);
515  break;
517  LinearAdvectionDiffusionReactionMatrixOp(inarray,outarray,mkey);
518  break;
519  case eLaplacian:
520  LaplacianMatrixOp(inarray,outarray,mkey);
521  break;
522  case eLaplacian00:
523  LaplacianMatrixOp(0,0,inarray,outarray,mkey);
524  break;
525  case eLaplacian01:
526  LaplacianMatrixOp(0,1,inarray,outarray,mkey);
527  break;
528  case eLaplacian02:
529  LaplacianMatrixOp(0,2,inarray,outarray,mkey);
530  break;
531  case eLaplacian10:
532  LaplacianMatrixOp(1,0,inarray,outarray,mkey);
533  break;
534  case eLaplacian11:
535  LaplacianMatrixOp(1,1,inarray,outarray,mkey);
536  break;
537  case eLaplacian12:
538  LaplacianMatrixOp(1,2,inarray,outarray,mkey);
539  break;
540  case eLaplacian20:
541  LaplacianMatrixOp(2,0,inarray,outarray,mkey);
542  break;
543  case eLaplacian21:
544  LaplacianMatrixOp(2,1,inarray,outarray,mkey);
545  break;
546  case eLaplacian22:
547  LaplacianMatrixOp(2,2,inarray,outarray,mkey);
548  break;
549  case eHelmholtz:
550  HelmholtzMatrixOp(inarray,outarray,mkey);
551  break;
552  default:
553  NEKERROR(ErrorUtil::efatal, "This matrix does not have an operator");
554  break;
555  }
556  }
557 
559  Array<OneD,NekDouble> &outarray,
560  const StdMatrixKey &mkey)
561  {
562  switch(mkey.GetMatrixType())
563  {
564  case eMass:
565  MassMatrixOp_MatFree(inarray,outarray,mkey);
566  break;
567  case eWeakDeriv0:
568  WeakDerivMatrixOp_MatFree(0,inarray,outarray,mkey);
569  break;
570  case eWeakDeriv1:
571  WeakDerivMatrixOp_MatFree(1,inarray,outarray,mkey);
572  break;
573  case eWeakDeriv2:
574  WeakDerivMatrixOp_MatFree(2,inarray,outarray,mkey);
575  break;
577  WeakDirectionalDerivMatrixOp_MatFree(inarray,outarray,mkey);
578  break;
579  case eMassLevelCurvature:
580  MassLevelCurvatureMatrixOp_MatFree(inarray,outarray,mkey);
581  break;
583  LinearAdvectionDiffusionReactionMatrixOp_MatFree(inarray,outarray,mkey,false);
584  break;
587  break;
588  case eLaplacian:
589  LaplacianMatrixOp_MatFree(inarray,outarray,mkey);
590  break;
591  case eLaplacian00:
592  LaplacianMatrixOp_MatFree(0,0,inarray,outarray,mkey);
593  break;
594  case eLaplacian01:
595  LaplacianMatrixOp_MatFree(0,1,inarray,outarray,mkey);
596  break;
597  case eLaplacian02:
598  LaplacianMatrixOp_MatFree(0,2,inarray,outarray,mkey);
599  break;
600  case eLaplacian10:
601  LaplacianMatrixOp_MatFree(1,0,inarray,outarray,mkey);
602  break;
603  case eLaplacian11:
604  LaplacianMatrixOp_MatFree(1,1,inarray,outarray,mkey);
605  break;
606  case eLaplacian12:
607  LaplacianMatrixOp_MatFree(1,2,inarray,outarray,mkey);
608  break;
609  case eLaplacian20:
610  LaplacianMatrixOp_MatFree(2,0,inarray,outarray,mkey);
611  break;
612  case eLaplacian21:
613  LaplacianMatrixOp_MatFree(2,1,inarray,outarray,mkey);
614  break;
615  case eLaplacian22:
616  LaplacianMatrixOp_MatFree(2,2,inarray,outarray,mkey);
617  break;
618  case eHelmholtz:
619  HelmholtzMatrixOp_MatFree(inarray,outarray,mkey);
620  break;
621  default:
622  NEKERROR(ErrorUtil::efatal, "This matrix does not have an operator");
623  break;
624  }
625  }
626 
628  Array<OneD,NekDouble> &outarray,
629  const StdMatrixKey &mkey)
630  {
631  int nq = GetTotPoints();
632  Array<OneD, NekDouble> tmp(nq);
633 
634  v_BwdTrans(inarray,tmp);
635 
636  if(mkey.HasVarCoeff(eVarCoeffMass))
637  {
638  Vmath::Vmul(nq, mkey.GetVarCoeff(eVarCoeffMass), 1, tmp, 1, tmp, 1);
639  }
640 
641  v_IProductWRTBase(tmp, outarray);
642  }
643 
644  void StdExpansion::LaplacianMatrixOp_MatFree(const int k1, const int k2,
645  const Array<OneD, const NekDouble> &inarray,
646  Array<OneD,NekDouble> &outarray,
647  const StdMatrixKey &mkey)
648  {
649  ASSERTL1(k1 >= 0 && k1 < GetCoordim(),"invalid first argument");
650  ASSERTL1(k2 >= 0 && k2 < GetCoordim(),"invalid second argument");
651 
652  int nq = GetTotPoints();
653  Array<OneD, NekDouble> tmp(nq);
654  Array<OneD, NekDouble> dtmp(nq);
655  VarCoeffType varcoefftypes[3][3]
659  };
660 
661  v_BwdTrans(inarray,tmp);
662  v_PhysDeriv(k2,tmp,dtmp);
663  if (mkey.GetNVarCoeff())
664  {
665  if (k1 == k2)
666  {
667  // By default, k1 == k2 has \sigma = 1 (diagonal entries)
668  if(mkey.HasVarCoeff(varcoefftypes[k1][k1]))
669  {
670  Vmath::Vmul(nq, mkey.GetVarCoeff(varcoefftypes[k1][k1]), 1, dtmp, 1, dtmp, 1);
671  }
672  v_IProductWRTDerivBase(k1, dtmp, outarray);
673  }
674  else
675  {
676  // By default, k1 != k2 has \sigma = 0 (off-diagonal entries)
677  if(mkey.HasVarCoeff(varcoefftypes[k1][k2]))
678  {
679  Vmath::Vmul(nq, mkey.GetVarCoeff(varcoefftypes[k1][k2]), 1, dtmp, 1, dtmp, 1);
680  v_IProductWRTDerivBase(k1, dtmp, outarray);
681  }
682  else
683  {
684  Vmath::Zero(GetNcoeffs(), outarray, 1);
685  }
686  }
687  }
688  else
689  {
690  // Multiply by svv tensor
692  {
693  Vmath::Vcopy(nq, dtmp, 1, tmp, 1);
694  SVVLaplacianFilter(dtmp,mkey);
695  Vmath::Vadd(nq, tmp, 1, dtmp, 1, dtmp, 1);
696  }
697  v_IProductWRTDerivBase(k1, dtmp, outarray);
698  }
699  }
700 
702  Array<OneD,NekDouble> &outarray,
703  const StdMatrixKey &mkey)
704  {
705  const int dim = GetCoordim();
706 
707  int i,j;
708 
710  Array<OneD,NekDouble> store2(m_ncoeffs,0.0);
711 
712  if(mkey.GetNVarCoeff() == 0)
713  {
714  // just call diagonal matrix form of laplcian operator
715  for(i = 0; i < dim; ++i)
716  {
717  LaplacianMatrixOp(i,i,inarray,store,mkey);
718  Vmath::Vadd(m_ncoeffs, store, 1, store2, 1, store2, 1);
719  }
720  }
721  else
722  {
723  const MatrixType mtype[3][3]
727  StdMatrixKeySharedPtr mkeyij;
728 
729  for(i = 0; i < dim; i++)
730  {
731  for(j = 0; j < dim; j++)
732  {
733  mkeyij = MemoryManager<StdMatrixKey>::AllocateSharedPtr(mkey,mtype[i][j]);
734  LaplacianMatrixOp(i,j,inarray,store,*mkeyij);
735  Vmath::Vadd(m_ncoeffs, store, 1, store2, 1, store2, 1);
736  }
737  }
738  }
739 
740  Vmath::Vcopy(m_ncoeffs,store2.get(),1,outarray.get(),1);
741  }
742 
744  const Array<OneD, const NekDouble> &inarray,
745  Array<OneD,NekDouble> &outarray,
746  const StdMatrixKey &mkey)
747  {
749  int nq = GetTotPoints();
750 
751  v_BwdTrans(inarray,tmp);
752  v_PhysDeriv(k1,tmp,tmp);
753 
755  if(mkey.HasVarCoeff(keys[k1]))
756  {
757  Vmath::Vmul(nq, &(mkey.GetVarCoeff(keys[k1]))[0], 1, &tmp[0], 1, &tmp[0], 1);
758  }
759 
760  v_IProductWRTBase(tmp, outarray);
761  }
762 
764  Array<OneD,NekDouble> &outarray,
765  const StdMatrixKey &mkey)
766  {
767  int nq = GetTotPoints();
768  // int varsize = ((mkey.GetVariableCoefficient(0)).num_elements())/dim;
769  Array<OneD, NekDouble> tmp(nq);
770 
771  v_BwdTrans(inarray,tmp);
772  // For Deformed mesh ==============
773  // if (varsize==nq)
774  // {
775  // v_PhysDirectionalDeriv(tmp,mkey.GetVariableCoefficient(0),tmp);
776  // }
777  //
778  // // For Regular mesh ==========
779  // else
780  // {
781  // ASSERTL0(false, "Wrong route");
782  // }
783 
784  v_IProductWRTBase(tmp, outarray);
785  }
786 
788  Array<OneD,NekDouble> &outarray,
789  const StdMatrixKey &mkey)
790  {
791  ///@todo fix this
792  // int nqtot = GetTotPoints();
793  // int matrixid = mkey.GetMatrixID();
794  //
795  // NekDouble checkweight=0.0;
796  // Array<OneD, NekDouble> tmp(nqtot), tan(nqtot), dtan0(nqtot), dtan1(nqtot), weight(nqtot,0.0);
797  //
798  // int gmatnumber = (mkey.GetVariableCoefficient(1)).num_elements();
799  //
800  // v_BwdTrans(inarray,tmp);
801  //
802  // // weight = \grad \cdot tanvec
803  // for(int k = 0; k < GetCoordim(); ++k)
804  // {
805  // Vmath::Vcopy(nqtot, &(mkey.GetVariableCoefficient(0))[k*nqtot],
806  // 1, &tan[0], 1);
807  //
808  // // For Regular mesh ...
809  // if(gmatnumber==1)
810  // {
811  // // D_{/xi} and D_{/eta}
812  // v_PhysDeriv(0,tan,dtan0);
813  // v_PhysDeriv(1,tan,dtan1);
814  //
815  // // d v / d x_i = (d \xi / d x_i)*( d v / d \xi ) + (d \eta / d x_i)*( d v / d \eta )
816  // Vmath::Svtvp(nqtot,(mkey.GetVariableCoefficient(2*k+1))[0],&dtan0[0],1,&weight[0],1,&weight[0],1);
817  // Vmath::Svtvp(nqtot,(mkey.GetVariableCoefficient(2*k+2))[0],&dtan1[0],1,&weight[0],1,&weight[0],1);
818  // }
819  //
820  // // For Curved mesh ...
821  // else if(gmatnumber==nqtot)
822  // {
823  // // D_{x} and D_{y}
824  // v_PhysDeriv(k,tan,dtan0);
825  // Vmath::Vadd(nqtot,&dtan0[0],1,&weight[0],1,&weight[0],1);
826  // }
827  //
828  // else
829  // {
830  // ASSERTL1( ((gmatnumber=1) || (gmatnumber==nqtot) ), "Gmat is not in a right size");
831  // }
832  // }
833  //
834  // Vmath::Vmul(nqtot, &weight[0], 1, &tmp[0], 1, &tmp[0], 1);
835  // v_IProductWRTBase(tmp, outarray);
836  }
837 
839  Array<OneD,NekDouble> &outarray,
840  const StdMatrixKey &mkey,
841  bool addDiffusionTerm)
842  {
843 
844  int i;
845  int ndir = mkey.GetNVarCoeff(); // assume num.r consts corresponds to directions
846  ASSERTL0(ndir,"Must define at least one advection velocity");
847 
848  NekDouble lambda = mkey.GetConstFactor(eFactorLambda);
849  int totpts = GetTotPoints();
850  Array<OneD, NekDouble> tmp(3*totpts);
851  Array<OneD, NekDouble> tmp_deriv = tmp + totpts;
852  Array<OneD, NekDouble> tmp_adv = tmp_deriv + totpts;
853 
854 
855  ASSERTL1(ndir <= GetCoordim(),"Number of constants is larger than coordinate dimensions");
856 
857  v_BwdTrans(inarray,tmp);
858 
859  VarCoeffType varcoefftypes[] = {eVarCoeffVelX, eVarCoeffVelY};
860 
861  //calculate u dx + v dy + ..
862  Vmath::Zero(totpts,tmp_adv,1);
863  for(i = 0; i < ndir; ++i)
864  {
865  v_PhysDeriv(i,tmp,tmp_deriv);
866  Vmath::Vvtvp(totpts,mkey.GetVarCoeff(varcoefftypes[i]),1,tmp_deriv,1,tmp_adv,1,tmp_adv,1);
867  }
868 
869  if(lambda) // add -lambda*u
870  {
871  Vmath::Svtvp(totpts,-lambda,tmp,1,tmp_adv,1,tmp_adv,1);
872  }
873 
874 
875  if(addDiffusionTerm)
876  {
878  StdMatrixKey mkeylap(eLaplacian,DetShapeType(),*this);
879  LaplacianMatrixOp(inarray,lap,mkeylap);
880 
881  v_IProductWRTBase(tmp_adv, outarray);
882  // Lap v - u.grad v + lambda*u
883  // => (grad u, grad v) + u.grad v - lambda*u
884  Vmath::Vadd(m_ncoeffs,lap,1,outarray,1,outarray,1);
885  }
886  else
887  {
888  v_IProductWRTBase(tmp_adv, outarray);
889  }
890 
891  }
892 
893 
895  Array<OneD,NekDouble> &outarray,
896  const StdMatrixKey &mkey)
897  {
898  NekDouble lambda = mkey.GetConstFactor(eFactorLambda);
900  StdMatrixKey mkeymass(eMass,DetShapeType(),*this);
901  StdMatrixKey mkeylap(eLaplacian,DetShapeType(),*this);
902 
903  MassMatrixOp(inarray,tmp,mkeymass);
904  LaplacianMatrixOp(inarray,outarray,mkeylap);
905 
906  Blas::Daxpy(m_ncoeffs, lambda, tmp, 1, outarray, 1);
907  }
908 
910  Array<OneD, NekDouble> &outarray)
911  {
912  int nq = GetTotPoints();
913  StdMatrixKey bwdtransmatkey(eBwdTrans,DetShapeType(),*this);
914  DNekMatSharedPtr bwdtransmat = GetStdMatrix(bwdtransmatkey);
915 
916  Blas::Dgemv('N',nq,m_ncoeffs,1.0,bwdtransmat->GetPtr().get(),
917  nq, inarray.get(), 1, 0.0, outarray.get(), 1);
918  }
919 
920  // VIRTUAL INLINE FUNCTIONS FROM HEADER FILE
921  void StdExpansion::SetUpPhysNormals(const int edge)
922  {
923  v_SetUpPhysNormals(edge);
924  }
925 
927  const Array<OneD, const NekDouble> &physvals)
928  {
929  return v_StdPhysEvaluate(Lcoord,physvals);
930  }
931 
933  {
934  return m_elmt_id;
935  }
936 
938  {
939  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
940  return NullNekDouble1DArray;
941  }
942 
943 
945  {
946  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
947  }
948 
949  void StdExpansion::v_SetUpPhysNormals(const int edge)
950  {
951  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
952  }
953 
954  int StdExpansion::v_CalcNumberOfCoefficients(const std::vector<unsigned int> &nummodes, int &modes_offset)
955  {
956  NEKERROR(ErrorUtil::efatal, "This function is not defined for this class");
957  return 0;
958  }
959 
961  const std::vector<unsigned int > &nummodes,
962  const int nmode_offset,
963  NekDouble *coeffs)
964  {
965  NEKERROR(ErrorUtil::efatal, "This function is not defined for this class");
966  }
967 
968 
969 
971  {
972  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
973  }
974 
979  Array< OneD, NekDouble> &outarray)
980  {
981  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
982  }
983 
985  {
986  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
988  }
989 
991  {
992  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
993  }
994 
996 
997  {
998  NEKERROR(ErrorUtil::efatal, "This function is only valid for three-dimensional LocalRegions");
1000  }
1001 
1003  {
1004  NEKERROR(ErrorUtil::efatal, "This function is only valid for two-dimensional LocalRegions");
1005  return eForwards;
1006  }
1007 
1009  {
1010  NEKERROR(ErrorUtil::efatal, "This function is only valid for one-dimensional LocalRegions");
1011  return eFwd;
1012  }
1013 
1014 
1016  {
1017  NEKERROR(ErrorUtil::efatal, "This function is only valid for two-dimensional LocalRegions");
1018  return eForwards;
1019  }
1020 
1021 
1024  Array<OneD, NekDouble> &outarray)
1025  {
1026  NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1027  }
1028 
1030  Array<OneD, NekDouble> &coeffs,
1032  {
1033  NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1034  }
1035 
1036 
1038  const Array<OneD, const NekDouble> &physvals)
1039 
1040  {
1041  NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1042  return 0;
1043  }
1044 
1046  {
1047  NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1048  }
1049 
1051  {
1052  ASSERTL0(false, "This function is needs defining for this shape");
1053  return 0;
1054  }
1055 
1057  {
1058  ASSERTL0(false, "This function is needs defining for this shape");
1059  return 0;
1060  }
1061 
1063  {
1064  ASSERTL0(false, "This function is needs defining for this shape");
1065  return 0;
1066  }
1067 
1069  {
1070  ASSERTL0(false, "This function is needs defining for this shape");
1071  return 0;
1072  }
1073 
1074  int StdExpansion::v_GetEdgeNcoeffs(const int i) const
1075  {
1076  ASSERTL0(false, "This function is not valid or not defined");
1077  return 0;
1078  }
1079 
1081  {
1082  ASSERTL0(false, "This function is not valid or not defined");
1083  return 0;
1084  }
1085 
1086  int StdExpansion::v_GetEdgeNumPoints(const int i) const
1087  {
1088  ASSERTL0(false, "This function is not valid or not defined");
1089  return 0;
1090  }
1091 
1093  {
1094  ASSERTL0(false, "This function is not valid or not defined");
1095  return 0;
1096  }
1097 
1099  {
1100  ASSERTL0(false, "This function is not valid or not defined");
1102  }
1103 
1104  const LibUtilities::BasisKey StdExpansion::v_DetFaceBasisKey(const int i, const int k) const
1105  {
1106  ASSERTL0(false, "This function is not valid or not defined");
1108  }
1109 
1110  int StdExpansion::v_GetFaceNumPoints(const int i) const
1111  {
1112  ASSERTL0(false, "This function is not valid or not defined");
1113  return 0;
1114  }
1115 
1116  int StdExpansion::v_GetFaceNcoeffs(const int i) const
1117  {
1118  ASSERTL0(false, "This function is not valid or not defined");
1119  return 0;
1120  }
1121 
1122  int StdExpansion::v_GetFaceIntNcoeffs(const int i) const
1123  {
1124  ASSERTL0(false, "This function is not valid or not defined");
1125  return 0;
1126  }
1127 
1129  {
1130  ASSERTL0(false, "This function is not valid or not defined");
1131  return 0;
1132  }
1133 
1135  {
1136  ASSERTL0(false, "This function is not valid or not defined");
1138  }
1139 
1141  {
1142  ASSERTL0(false, "This function is not valid or not defined");
1143 
1145  }
1146 
1148  {
1149  ASSERTL0(false, "This function is not valid or not defined");
1150 
1152  }
1153 
1155  {
1156  ASSERTL0(false, "This expansion does not have a shape type defined");
1158  }
1159 
1160  boost::shared_ptr<StdExpansion>
1162  {
1163  ASSERTL0(false,"This method is not defined for this expansion");
1164  StdExpansionSharedPtr returnval;
1165  return returnval;
1166  }
1167 
1169  {
1170  ASSERTL0(false, "This function is not valid or not defined");
1171  return 0;
1172  }
1173 
1175  {
1176  ASSERTL0(false,"This function has not been defined for this expansion");
1177  return false;
1178  }
1179 
1180 
1182  {
1183  return false;
1184  }
1185 
1187  const Array<OneD, const NekDouble>& inarray,
1188  Array<OneD, NekDouble> &outarray)
1189  {
1190  NEKERROR(ErrorUtil::efatal, "This method has not been defined");
1191  }
1192 
1193  /**
1194  *
1195  */
1197  Array<OneD, NekDouble> &outarray)
1198  {
1199  NEKERROR(ErrorUtil::efatal, "This method has not been defined");
1200  }
1201 
1202 
1203  /**
1204  * @brief Integrates the specified function over the domain.
1205  * @see StdRegions#StdExpansion#Integral.
1206  */
1208  {
1209  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1210  "local expansions");
1211  return 0;
1212  }
1213 
1214 
1215  void StdExpansion::v_AddRobinMassMatrix(const int edgeid, const Array<OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
1216  {
1217  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1218  "specific element types");
1219  }
1220 
1222  {
1223  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1224  "specific element types");
1225  }
1226 
1227  /**
1228  * @brief Calculate the derivative of the physical points
1229  * @see StdRegions#StdExpansion#PhysDeriv
1230  */
1232  Array<OneD, NekDouble> &out_d1,
1233  Array<OneD, NekDouble> &out_d2,
1234  Array<OneD, NekDouble> &out_d3)
1235  {
1236  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1237  "local expansions");
1238  }
1239 
1241  Array<OneD, NekDouble> &out_ds)
1242  {
1243  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1244  "local expansions");
1245  }
1247  Array<OneD, NekDouble>& out_dn)
1248  {
1249  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1250  "local expansions");
1251  }
1252 
1253  /**
1254  * @brief Calculate the derivative of the physical points in a
1255  * given direction
1256  * @see StdRegions#StdExpansion#PhysDeriv
1257  */
1258  void StdExpansion::v_PhysDeriv(const int dir,
1259  const Array<OneD, const NekDouble>& inarray,
1260  Array<OneD, NekDouble> &out_d0)
1261 
1262  {
1263  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1264  "specific element types");
1265  }
1266 
1267  /**
1268  * @brief Physical derivative along a direction vector.
1269  * @see StdRegions#StdExpansion#PhysDirectionalDeriv
1270  */
1272  const Array<OneD, const NekDouble>& direction,
1273  Array<OneD, NekDouble> &outarray)
1274  {
1275  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1276  "specific element types");
1277  }
1278 
1280  Array<OneD, NekDouble> &out_d1,
1281  Array<OneD, NekDouble> &out_d2,
1282  Array<OneD, NekDouble> &out_d3)
1283  {
1284  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1285  }
1286 
1287  void StdExpansion::v_StdPhysDeriv (const int dir,
1288  const Array<OneD, const NekDouble>& inarray,
1289  Array<OneD, NekDouble> &outarray)
1290  {
1291  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1292  }
1293 
1294 
1296  {
1297  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1298  return 0;
1299  }
1300 
1301 
1303  {
1304  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1305  return 0;
1306  }
1307 
1308 
1309  void StdExpansion::v_FillMode(const int mode, Array<OneD, NekDouble> &outarray)
1310  {
1311  NEKERROR(ErrorUtil::efatal, "This function has not "
1312  "been defined for this shape");
1313  }
1314 
1316  {
1317  NEKERROR(ErrorUtil::efatal, "This function has not "
1318  "been defined for this element");
1319  DNekMatSharedPtr returnval;
1320  return returnval;
1321  }
1322 
1324  {
1325  NEKERROR(ErrorUtil::efatal, "This function has not "
1326  "been defined for this element");
1327  DNekMatSharedPtr returnval;
1328  return returnval;
1329  }
1330 
1332  Array<OneD, NekDouble> &coords_1,
1333  Array<OneD, NekDouble> &coords_2)
1334  {
1335  NEKERROR(ErrorUtil::efatal, "Write coordinate definition method");
1336  }
1337 
1339  Array<OneD, NekDouble> &coord)
1340  {
1341  NEKERROR(ErrorUtil::efatal, "Write coordinate definition method");
1342  }
1343 
1345  {
1346  NEKERROR(ErrorUtil::efatal, "Write method");
1347  return -1;
1348  }
1349 
1351  {
1352  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1353  }
1354 
1356  {
1357  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1358  }
1359 
1360  int StdExpansion::v_GetVertexMap(const int localVertexId,
1361  bool useCoeffPacking)
1362  {
1363  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1364  return 0;
1365  }
1366 
1367  void StdExpansion::v_GetEdgeInteriorMap(const int eid, const Orientation edgeOrient,
1368  Array<OneD, unsigned int> &maparray,
1369  Array<OneD, int> &signarray)
1370  {
1371  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1372  }
1373 
1374  void StdExpansion::v_GetFaceInteriorMap(const int fid, const Orientation faceOrient,
1375  Array<OneD, unsigned int> &maparray,
1376  Array<OneD, int> &signarray)
1377  {
1378  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1379  }
1380 
1381  void StdExpansion::v_GetEdgeToElementMap(const int eid, const Orientation edgeOrient,
1382  Array<OneD, unsigned int> &maparray,
1383  Array<OneD, int> &signarray)
1384  {
1385  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1386  }
1387 
1388  void StdExpansion::v_GetFaceToElementMap(const int fid, const Orientation faceOrient,
1389  Array<OneD, unsigned int> &maparray,
1390  Array<OneD, int> &signarray,
1391  int nummodesA, int nummodesB)
1392  {
1393  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1394  }
1395 
1397  {
1398  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape or library" );
1399  }
1400 
1401  void StdExpansion::v_GetEdgePhysVals(const int edge, const boost::shared_ptr<StdExpansion> &EdgeExp, const Array<OneD, const NekDouble> &inarray, Array<OneD,NekDouble> &outarray)
1402  {
1403  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape or library" );
1404  }
1405 
1406  void StdExpansion::v_GetTracePhysVals(const int edge, const boost::shared_ptr<StdExpansion> &EdgeExp, const Array<OneD, const NekDouble> &inarray, Array<OneD,NekDouble> &outarray, StdRegions::Orientation orient)
1407  {
1408  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape or library" );
1409  }
1410 
1411  void StdExpansion::v_GetVertexPhysVals(const int vertex, const Array<OneD, const NekDouble> &inarray, NekDouble &outarray)
1412  {
1413  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape or library" );
1414  }
1415 
1417  {
1418  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape or library" );
1419  }
1420 
1422  const int edge,
1423  Array<OneD, NekDouble> &outarray)
1424  {
1426  "Method does not exist for this shape or library");
1427  }
1428 
1429  void StdExpansion::v_GetFacePhysVals( const int face,
1430  const boost::shared_ptr<StdExpansion> &FaceExp,
1431  const Array<OneD, const NekDouble> &inarray,
1432  Array<OneD, NekDouble> &outarray,
1433  StdRegions::Orientation orient)
1434  {
1435  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape or library" );
1436  }
1437 
1439  const Array<OneD, const NekDouble> &inarray,
1440  Array<OneD, NekDouble> &outarray)
1441  {
1442  v_MultiplyByStdQuadratureMetric(inarray,outarray);
1443  }
1444 
1446  const Array<OneD, const NekDouble> &inarray,
1447  Array<OneD, NekDouble> &outarray)
1448  {
1449  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape or library");
1450  }
1451 
1452  const boost::shared_ptr<SpatialDomains::GeomFactors>& StdExpansion::v_GetMetricInfo() const
1453  {
1454  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
1456 
1457  }
1458 
1460  Array<OneD, NekDouble> &outarray)
1461  {
1462  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1463  }
1464 
1466  Array<OneD, NekDouble> &outarray,
1467  bool multiplybyweights)
1468  {
1469  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1470  }
1471 
1473  const Array<OneD, const NekDouble>& inarray,
1474  Array<OneD, NekDouble> &outarray)
1475  {
1476  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1477  }
1478 
1480  Array<OneD,NekDouble> &outarray,
1481  const StdMatrixKey &mkey)
1482  {
1483  // If this function is not reimplemented on shape level, the function
1484  // below will be called
1485  MassMatrixOp_MatFree(inarray,outarray,mkey);
1486  }
1487 
1489  Array<OneD,NekDouble> &outarray,
1490  const StdMatrixKey &mkey)
1491  {
1492  // If this function is not reimplemented on shape level, the function
1493  // below will be called
1494  LaplacianMatrixOp_MatFree(inarray,outarray,mkey);
1495  }
1496 
1497 
1499  const StdMatrixKey &mkey)
1500  {
1501  ASSERTL0(false, "This function is not defined in StdExpansion.");
1502  }
1503 
1505  const Array<OneD, const NekDouble> &inarray,
1506  Array<OneD, NekDouble> &outarray)
1507  {
1508  ASSERTL0(false, "This function is not defined in StdExpansion.");
1509  }
1510 
1511  void StdExpansion::v_LaplacianMatrixOp(const int k1, const int k2,
1512  const Array<OneD, const NekDouble> &inarray,
1513  Array<OneD,NekDouble> &outarray,
1514  const StdMatrixKey &mkey)
1515  {
1516  // If this function is not reimplemented on shape level, the function
1517  // below will be called
1518  LaplacianMatrixOp_MatFree(k1,k2,inarray,outarray,mkey);
1519  }
1520 
1522  const Array<OneD, const NekDouble> &inarray,
1523  Array<OneD,NekDouble> &outarray,
1524  const StdMatrixKey &mkey)
1525  {
1526  // If this function is not reimplemented on shape level, the function
1527  // below will be called
1528  WeakDerivMatrixOp_MatFree(i,inarray,outarray,mkey);
1529 
1530  }
1531 
1533  Array<OneD,NekDouble> &outarray,
1534  const StdMatrixKey &mkey)
1535  {
1536  // If this function is not reimplemented on shape level, the function
1537  // below will be called
1538  WeakDirectionalDerivMatrixOp_MatFree(inarray,outarray,mkey);
1539 
1540  }
1541 
1543  Array<OneD,NekDouble> &outarray,
1544  const StdMatrixKey &mkey)
1545  {
1546  // If this function is not reimplemented on shape level, the function
1547  // below will be called
1548  MassLevelCurvatureMatrixOp_MatFree(inarray,outarray,mkey);
1549  }
1550 
1552  const NekDouble> &inarray,
1553  Array<OneD,NekDouble> &outarray,
1554  const StdMatrixKey &mkey, bool addDiffusionTerm)
1555  {
1556  // If this function is not reimplemented on shape level, the function
1557  // below will be called
1558  LinearAdvectionDiffusionReactionMatrixOp_MatFree(inarray,outarray,mkey,addDiffusionTerm);
1559 
1560  }
1561 
1563  Array<OneD,NekDouble> &outarray,
1564  const StdMatrixKey &mkey)
1565  {
1566  // If this function is not reimplemented on shape level, the function
1567  // below will be called
1568  HelmholtzMatrixOp_MatFree(inarray,outarray,mkey);
1569  }
1570 
1572  Array<OneD,NekDouble> &outarray,
1573  const StdMatrixKey &mkey)
1574  {
1575  // If this function is not reimplemented on shape level, the function
1576  // below will be called
1577  LaplacianMatrixOp_MatFree_GenericImpl(inarray,outarray,mkey);
1578  }
1579 
1581  const Array<OneD, const NekDouble> &inarray,
1582  Array<OneD, NekDouble> &outarray,
1584  {
1585  ASSERTL0(false, "Not implemented.");
1586  }
1587 
1589  Array<OneD,NekDouble> &outarray,
1590  const StdMatrixKey &mkey)
1591  {
1592  // If this function is not reimplemented on shape level, the function
1593  // below will be called
1594  HelmholtzMatrixOp_MatFree_GenericImpl(inarray,outarray,mkey);
1595  }
1596 
1597  const NormalVector & StdExpansion::v_GetEdgeNormal(const int edge) const
1598  {
1599  ASSERTL0(false, "Cannot get edge normals for this expansion.");
1600  static NormalVector result;
1601  return result;
1602  }
1603 
1605  {
1606  ASSERTL0(false, "Cannot compute edge normal for this expansion.");
1607  }
1608 
1610  {
1611  ASSERTL0(false, "Not implemented.");
1612  }
1613 
1615  {
1616  ASSERTL0(false, "Not implemented.");
1617  return false;
1618  }
1619 
1621  {
1622  ASSERTL0(false, "Cannot compute face normal for this expansion.");
1623  }
1624 
1626  {
1627  ASSERTL0(false, "Not implemented.");
1628  }
1629 
1631  {
1632  ASSERTL0(false, "Cannot compute vertex normal for this expansion.");
1633  }
1634 
1635  const NormalVector & StdExpansion::v_GetFaceNormal(const int face) const
1636  {
1637  ASSERTL0(false, "Cannot get face normals for this expansion.");
1638  static NormalVector result;
1639  return result;
1640  }
1641 
1642  const NormalVector & StdExpansion::v_GetVertexNormal(const int vertex) const
1643  {
1644  ASSERTL0(false, "Cannot get vertex normals for this expansion.");
1645  static NormalVector result;
1646  return result;
1647  }
1648 
1650  {
1651  ASSERTL0(false, "Cannot get face normals for this expansion.");
1652  static NormalVector result;
1653  return result;
1654  }
1655 
1658  {
1659  ASSERTL0(false, "Not implemented.");
1660  Array<OneD, unsigned int> noinversemap(1);
1661  return noinversemap;
1662  }
1663 
1666  StdRegions::Orientation faceOrient)
1667  {
1668  ASSERTL0(false, "Not implemented.");
1669  Array<OneD, unsigned int> noinversemap(1);
1670  return noinversemap;
1671  }
1672 
1675  const DNekScalMatSharedPtr & m_transformationmatrix)
1676  {
1677  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
1678  return NullDNekMatSharedPtr;
1679  }
1680 
1682  const Array<OneD, const NekDouble> &inarray,
1683  Array<OneD, NekDouble> &outarray)
1684  {
1686  StdMatrixKey Ikey(ePhysInterpToEquiSpaced, shape, *this);
1687  DNekMatSharedPtr intmat = GetStdMatrix(Ikey);
1688 
1689  int nqtot = 1;
1690  int nqbase;
1691  int np = 0;
1692  for(int i = 0; i < m_base.num_elements(); ++i)
1693  {
1694  nqbase = m_base[i]->GetNumPoints();
1695  nqtot *= nqbase;
1696  np = max(np,nqbase);
1697  }
1698 
1699  NekVector<NekDouble> in (nqtot,inarray,eWrapper);
1701  out = (*intmat) * in;
1702  }
1703 
1705  Array<OneD, int> &conn,
1706  bool standard)
1707  {
1708  ASSERTL0(false, "Not implemented.");
1709  }
1710  }//end namespace
1711 }//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)
static const PointsKey NullPointsKey(0, eNoPointsType)
virtual int v_GetVertexMap(int localVertexId, bool useCoeffPacking=false)
virtual void v_ComputeFaceNormal(const int face)
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:454
NekDouble GetConstFactor(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:122
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
virtual void v_NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
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:158
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 int v_GetTotalEdgeIntNcoeffs() const
boost::shared_ptr< IndexMapValues > IndexMapValuesSharedPtr
Definition: IndexMapKey.h:126
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:82
virtual int v_NumBndryCoeffs() const
static boost::shared_ptr< GeomFactors > NullGeomFactorsSharedPtr
virtual void v_NegateFaceNormal(const int face)
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:931
static Array< OneD, NekDouble > NullNekDouble1DArray
virtual void v_ExtractDataToCoeffs(const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmode_offset, NekDouble *coeffs)
Unpack data from input file assuming it comes from the same expansion type.
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
virtual void v_FwdTrans_BndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
bool HasVarCoeff(const StdRegions::VarCoeffType &coeff) const
Definition: StdMatrixKey.h:173
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:22
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)
virtual LibUtilities::ShapeType v_DetShapeType() const
virtual const NormalVector & v_GetVertexNormal(const int vertex) const
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)
virtual void v_GetFaceInteriorMap(const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
virtual int v_GetShapeDimension() const
void WeakDirectionalDerivMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual const NormalVector & v_GetEdgeNormal(const int edge) const
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:952
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)
const Array< OneD, const NekDouble > & GetVarCoeff(const StdRegions::VarCoeffType &coeff) const
Definition: StdMatrixKey.h:152
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:471
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:428
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
STL namespace.
virtual int v_NumDGBndryCoeffs() const
static const NekDouble kNekSqrtTol
virtual void v_DropLocStaticCondMatrix(const LocalRegions::MatrixKey &mkey)
virtual ~StdExpansion()
Destructor.
LibUtilities::ShapeType GetShapeType() const
Definition: StdMatrixKey.h:87
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.
virtual int v_GetEdgeNcoeffs(const int i) const
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)
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
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.
void Vabs(int n, const T *x, const int incx, T *y, const int incy)
vabs: y = |x|
Definition: Vmath.cpp:410
virtual const LibUtilities::PointsKey v_GetNodalPointsKey() const
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:684
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
void IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
void PhysInterpToSimplexEquiSpaced(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs an interpolation from the physical space points provided at input into an arra...
int GetNumberOfCoefficients(ShapeType shape, std::vector< unsigned int > &modes, int offset)
Definition: ShapeType.hpp:312
static DNekMatSharedPtr NullDNekMatSharedPtr
Definition: NekTypeDefs.hpp:79
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_AddRobinMassMatrix(const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
virtual int v_GetFaceNcoeffs(const int i) const
void WeakDerivMatrixOp(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:966
void MassLevelCurvatureMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:981
bool ConstFactorExists(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:131
BasisManagerT & BasisManager(void)
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:141
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.
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)
virtual LibUtilities::BasisType v_GetEdgeBasisType(const int i) const
The base class for all shapes.
Definition: StdExpansion.h:69
virtual int v_GetEdgeNumPoints(const int i) const
virtual void SetUpPhysNormals(const int edge)
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 const NormalVector & v_GetSurfaceNormal(const int id) const
LibUtilities::PointsType GetNodalPointsType() const
Definition: StdMatrixKey.h:92
virtual void v_GetCoord(const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey)
IndexMapValuesSharedPtr CreateIndexMap(const IndexMapKey &ikey)
Create an IndexMap which contains mapping information linking any specific element shape with either ...
virtual LibUtilities::PointsKey v_GetFacePointsKey(const int i, const int j) const
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 LibUtilities::BasisKey v_DetEdgeBasisKey(const int i) 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)
boost::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:74
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)
virtual boost::shared_ptr< StdExpansion > v_GetStdExp(void) const
virtual int v_GetFaceIntNcoeffs(const int i) const
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:795
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
virtual int v_GetTotalFaceIntNcoeffs() const
virtual void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
Defines a specification for a set of points.
Definition: Points.h:58
virtual const LibUtilities::BasisKey v_DetFaceBasisKey(const int i, const int k) const
double NekDouble
virtual bool v_IsBoundaryInteriorExpansion()
void BwdTrans_MatOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_GetTracePhysVals(const int edge, const boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=eNoOrientation)
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:802
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 Array< OneD, unsigned int > v_GetFaceInverseBoundaryMap(int fid, StdRegions::Orientation faceOrient=eNoOrientation)
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)
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.
void MassLevelCurvatureMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual const NormalVector & v_GetFaceNormal(const int face) const
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:329
boost::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:72
virtual StdRegions::Orientation v_GetCartesianEorient(int edge)
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:938
virtual void v_SetCoeffsToOrientation(StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual const boost::shared_ptr< SpatialDomains::GeomFactors > & v_GetMetricInfo() const
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual int v_GetNedges() const
virtual void v_ComputeVertexNormal(const int vertex)
virtual NekDouble v_StdPhysEvaluate(const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:213
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)
Orientation GetIndexOrientation() const
Definition: IndexMapKey.h:93
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:81
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)
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:131
virtual int v_GetFaceNumPoints(const int i) const
StdExpansion()
Default Constructor.
IndexMapType GetIndexMapType() const
Definition: IndexMapKey.h:88
void WeakDirectionalDerivMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:974
virtual void v_GetEdgeToElementMap(const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
virtual int v_GetNfaces() const
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:359
boost::shared_ptr< StdExpansion > StdExpansionSharedPtr
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void GeneralMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1038
virtual void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
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 StdRegions::Orientation v_GetPorient(int point)
virtual void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void GetBoundaryMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:790
boost::shared_ptr< StdMatrixKey > StdMatrixKeySharedPtr
Definition: StdMatrixKey.h:196
virtual DNekMatSharedPtr v_BuildInverseTransformationMatrix(const DNekScalMatSharedPtr &m_transformationmatrix)
Describes the specification for a Basis.
Definition: Basis.h:50
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:285
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:169
void BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:226
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:249
void LinearAdvectionDiffusionReactionMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
Definition: StdExpansion.h:988
virtual void v_GetFacePhysVals(const int face, const boost::shared_ptr< StdExpansion > &FaceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient)