Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 eEquiSpacedToCoeffs:
447  {
448  // check to see if equispaced basis
449  int nummodes = m_base[0]->GetNumModes();
450  bool equispaced = true;
451  for(int i = 1; i < m_base.num_elements(); ++i)
452  {
453  if(m_base[i]->GetNumModes() != nummodes)
454  {
455  equispaced = false;
456  }
457  }
458 
459  ASSERTL0(equispaced,
460  "Currently need to have same num modes in all "
461  "directionmodes to use EquiSpacedToCoeff method");
462 
463  int ntot = GetTotPoints();
464  Array<OneD, NekDouble> qmode(ntot);
466 
469  int cnt = 0;
470  for(int i = 0; i < m_ncoeffs; ++i)
471  {
472  // Get mode at quadrature points
473  FillMode(i,qmode);
474 
475  // interpolate to equi spaced
476  PhysInterpToSimplexEquiSpaced(qmode,emode,nummodes);
477 
478  // fill matrix
479  Vmath::Vcopy(m_ncoeffs, &emode[0], 1,
480  returnval->GetRawPtr() + i*m_ncoeffs, 1);
481  }
482  // invert matrix
483  returnval->Invert();
484 
485  }
486  break;
487  case eMass:
488  case eHelmholtz:
489  case eLaplacian:
490  case eLaplacian00:
491  case eLaplacian01:
492  case eLaplacian02:
493  case eLaplacian11:
494  case eLaplacian12:
495  case eLaplacian22:
496  case eWeakDeriv0:
497  case eWeakDeriv1:
498  case eWeakDeriv2:
500  case eMassLevelCurvature:
503  {
506  DNekMat &Mat = *returnval;
507 
508  for(i=0; i < m_ncoeffs; ++i)
509  {
510  Vmath::Zero(m_ncoeffs, tmp, 1);
511  tmp[i] = 1.0;
512 
513  GeneralMatrixOp_MatFree(tmp,tmp,mkey);
514 
515  Vmath::Vcopy(m_ncoeffs,&tmp[0],1,
516  &(Mat.GetPtr())[0]+i*m_ncoeffs,1);
517  }
518  }
519  break;
520  default:
521  {
522  NEKERROR(ErrorUtil::efatal, "This type of matrix can not be created using a general approach");
523  }
524  break;
525  }
526 
527  return returnval;
528  }
529 
531  Array<OneD,NekDouble> &outarray,
532  const StdMatrixKey &mkey)
533  {
534  switch(mkey.GetMatrixType())
535  {
536  case eMass:
537  MassMatrixOp(inarray,outarray,mkey);
538  break;
539  case eWeakDeriv0:
540  WeakDerivMatrixOp(0,inarray,outarray,mkey);
541  break;
542  case eWeakDeriv1:
543  WeakDerivMatrixOp(1,inarray,outarray,mkey);
544  break;
545  case eWeakDeriv2:
546  WeakDerivMatrixOp(2,inarray,outarray,mkey);
547  break;
549  WeakDirectionalDerivMatrixOp(inarray,outarray,mkey);
550  break;
551  case eMassLevelCurvature:
552  MassLevelCurvatureMatrixOp(inarray,outarray,mkey);
553  break;
555  LinearAdvectionDiffusionReactionMatrixOp(inarray,outarray,mkey,false);
556  break;
558  LinearAdvectionDiffusionReactionMatrixOp(inarray,outarray,mkey);
559  break;
560  case eLaplacian:
561  LaplacianMatrixOp(inarray,outarray,mkey);
562  break;
563  case eLaplacian00:
564  LaplacianMatrixOp(0,0,inarray,outarray,mkey);
565  break;
566  case eLaplacian01:
567  LaplacianMatrixOp(0,1,inarray,outarray,mkey);
568  break;
569  case eLaplacian02:
570  LaplacianMatrixOp(0,2,inarray,outarray,mkey);
571  break;
572  case eLaplacian10:
573  LaplacianMatrixOp(1,0,inarray,outarray,mkey);
574  break;
575  case eLaplacian11:
576  LaplacianMatrixOp(1,1,inarray,outarray,mkey);
577  break;
578  case eLaplacian12:
579  LaplacianMatrixOp(1,2,inarray,outarray,mkey);
580  break;
581  case eLaplacian20:
582  LaplacianMatrixOp(2,0,inarray,outarray,mkey);
583  break;
584  case eLaplacian21:
585  LaplacianMatrixOp(2,1,inarray,outarray,mkey);
586  break;
587  case eLaplacian22:
588  LaplacianMatrixOp(2,2,inarray,outarray,mkey);
589  break;
590  case eHelmholtz:
591  HelmholtzMatrixOp(inarray,outarray,mkey);
592  break;
593  default:
594  NEKERROR(ErrorUtil::efatal, "This matrix does not have an operator");
595  break;
596  }
597  }
598 
600  Array<OneD,NekDouble> &outarray,
601  const StdMatrixKey &mkey)
602  {
603  switch(mkey.GetMatrixType())
604  {
605  case eMass:
606  MassMatrixOp_MatFree(inarray,outarray,mkey);
607  break;
608  case eWeakDeriv0:
609  WeakDerivMatrixOp_MatFree(0,inarray,outarray,mkey);
610  break;
611  case eWeakDeriv1:
612  WeakDerivMatrixOp_MatFree(1,inarray,outarray,mkey);
613  break;
614  case eWeakDeriv2:
615  WeakDerivMatrixOp_MatFree(2,inarray,outarray,mkey);
616  break;
618  WeakDirectionalDerivMatrixOp_MatFree(inarray,outarray,mkey);
619  break;
620  case eMassLevelCurvature:
621  MassLevelCurvatureMatrixOp_MatFree(inarray,outarray,mkey);
622  break;
624  LinearAdvectionDiffusionReactionMatrixOp_MatFree(inarray,outarray,mkey,false);
625  break;
628  break;
629  case eLaplacian:
630  LaplacianMatrixOp_MatFree(inarray,outarray,mkey);
631  break;
632  case eLaplacian00:
633  LaplacianMatrixOp_MatFree(0,0,inarray,outarray,mkey);
634  break;
635  case eLaplacian01:
636  LaplacianMatrixOp_MatFree(0,1,inarray,outarray,mkey);
637  break;
638  case eLaplacian02:
639  LaplacianMatrixOp_MatFree(0,2,inarray,outarray,mkey);
640  break;
641  case eLaplacian10:
642  LaplacianMatrixOp_MatFree(1,0,inarray,outarray,mkey);
643  break;
644  case eLaplacian11:
645  LaplacianMatrixOp_MatFree(1,1,inarray,outarray,mkey);
646  break;
647  case eLaplacian12:
648  LaplacianMatrixOp_MatFree(1,2,inarray,outarray,mkey);
649  break;
650  case eLaplacian20:
651  LaplacianMatrixOp_MatFree(2,0,inarray,outarray,mkey);
652  break;
653  case eLaplacian21:
654  LaplacianMatrixOp_MatFree(2,1,inarray,outarray,mkey);
655  break;
656  case eLaplacian22:
657  LaplacianMatrixOp_MatFree(2,2,inarray,outarray,mkey);
658  break;
659  case eHelmholtz:
660  HelmholtzMatrixOp_MatFree(inarray,outarray,mkey);
661  break;
662  default:
663  NEKERROR(ErrorUtil::efatal, "This matrix does not have an operator");
664  break;
665  }
666  }
667 
669  Array<OneD,NekDouble> &outarray,
670  const StdMatrixKey &mkey)
671  {
672  int nq = GetTotPoints();
673  Array<OneD, NekDouble> tmp(nq);
674 
675  v_BwdTrans(inarray,tmp);
676 
677  if(mkey.HasVarCoeff(eVarCoeffMass))
678  {
679  Vmath::Vmul(nq, mkey.GetVarCoeff(eVarCoeffMass), 1, tmp, 1, tmp, 1);
680  }
681 
682  v_IProductWRTBase(tmp, outarray);
683  }
684 
685  void StdExpansion::LaplacianMatrixOp_MatFree(const int k1, const int k2,
686  const Array<OneD, const NekDouble> &inarray,
687  Array<OneD,NekDouble> &outarray,
688  const StdMatrixKey &mkey)
689  {
690  ASSERTL1(k1 >= 0 && k1 < GetCoordim(),"invalid first argument");
691  ASSERTL1(k2 >= 0 && k2 < GetCoordim(),"invalid second argument");
692 
693  int nq = GetTotPoints();
694  Array<OneD, NekDouble> tmp(nq);
695  Array<OneD, NekDouble> dtmp(nq);
696  VarCoeffType varcoefftypes[3][3]
700  };
701 
702  v_BwdTrans(inarray,tmp);
703  v_PhysDeriv(k2,tmp,dtmp);
704  if (mkey.GetNVarCoeff()&&
706  {
707  if (k1 == k2)
708  {
709  // By default, k1 == k2 has \sigma = 1 (diagonal entries)
710  if(mkey.HasVarCoeff(varcoefftypes[k1][k1]))
711  {
712  Vmath::Vmul(nq, mkey.GetVarCoeff(varcoefftypes[k1][k1]), 1, dtmp, 1, dtmp, 1);
713  }
714  v_IProductWRTDerivBase(k1, dtmp, outarray);
715  }
716  else
717  {
718  // By default, k1 != k2 has \sigma = 0 (off-diagonal entries)
719  if(mkey.HasVarCoeff(varcoefftypes[k1][k2]))
720  {
721  Vmath::Vmul(nq, mkey.GetVarCoeff(varcoefftypes[k1][k2]), 1, dtmp, 1, dtmp, 1);
722  v_IProductWRTDerivBase(k1, dtmp, outarray);
723  }
724  else
725  {
726  Vmath::Zero(GetNcoeffs(), outarray, 1);
727  }
728  }
729  }
730  else
731  {
732  // Multiply by svv tensor
734  {
735  Vmath::Vcopy(nq, dtmp, 1, tmp, 1);
736  SVVLaplacianFilter(dtmp,mkey);
737  Vmath::Vadd(nq, tmp, 1, dtmp, 1, dtmp, 1);
738  }
739  v_IProductWRTDerivBase(k1, dtmp, outarray);
740  }
741  }
742 
744  Array<OneD,NekDouble> &outarray,
745  const StdMatrixKey &mkey)
746  {
747  const int dim = GetCoordim();
748 
749  int i,j;
750 
752  Array<OneD,NekDouble> store2(m_ncoeffs,0.0);
753 
754  if(mkey.GetNVarCoeff() == 0||mkey.ConstFactorExists(eFactorSVVDiffCoeff))
755  {
756  // just call diagonal matrix form of laplcian operator
757  for(i = 0; i < dim; ++i)
758  {
759  LaplacianMatrixOp(i,i,inarray,store,mkey);
760  Vmath::Vadd(m_ncoeffs, store, 1, store2, 1, store2, 1);
761  }
762  }
763  else
764  {
765  const MatrixType mtype[3][3]
769  StdMatrixKeySharedPtr mkeyij;
770 
771  for(i = 0; i < dim; i++)
772  {
773  for(j = 0; j < dim; j++)
774  {
775  mkeyij = MemoryManager<StdMatrixKey>::AllocateSharedPtr(mkey,mtype[i][j]);
776  LaplacianMatrixOp(i,j,inarray,store,*mkeyij);
777  Vmath::Vadd(m_ncoeffs, store, 1, store2, 1, store2, 1);
778  }
779  }
780  }
781 
782  Vmath::Vcopy(m_ncoeffs,store2.get(),1,outarray.get(),1);
783  }
784 
786  const Array<OneD, const NekDouble> &inarray,
787  Array<OneD,NekDouble> &outarray,
788  const StdMatrixKey &mkey)
789  {
791  int nq = GetTotPoints();
792 
793  v_BwdTrans(inarray,tmp);
794  v_PhysDeriv(k1,tmp,tmp);
795 
797  if(mkey.HasVarCoeff(keys[k1]))
798  {
799  Vmath::Vmul(nq, &(mkey.GetVarCoeff(keys[k1]))[0], 1, &tmp[0], 1, &tmp[0], 1);
800  }
801 
802  v_IProductWRTBase(tmp, outarray);
803  }
804 
806  Array<OneD,NekDouble> &outarray,
807  const StdMatrixKey &mkey)
808  {
809  int nq = GetTotPoints();
810  // int varsize = ((mkey.GetVariableCoefficient(0)).num_elements())/dim;
811  Array<OneD, NekDouble> tmp(nq);
812 
813  v_BwdTrans(inarray,tmp);
814  // For Deformed mesh ==============
815  // if (varsize==nq)
816  // {
817  // v_PhysDirectionalDeriv(tmp,mkey.GetVariableCoefficient(0),tmp);
818  // }
819  //
820  // // For Regular mesh ==========
821  // else
822  // {
823  // ASSERTL0(false, "Wrong route");
824  // }
825 
826  v_IProductWRTBase(tmp, outarray);
827  }
828 
830  Array<OneD,NekDouble> &outarray,
831  const StdMatrixKey &mkey)
832  {
833  ///@todo fix this
834  // int nqtot = GetTotPoints();
835  // int matrixid = mkey.GetMatrixID();
836  //
837  // NekDouble checkweight=0.0;
838  // Array<OneD, NekDouble> tmp(nqtot), tan(nqtot), dtan0(nqtot), dtan1(nqtot), weight(nqtot,0.0);
839  //
840  // int gmatnumber = (mkey.GetVariableCoefficient(1)).num_elements();
841  //
842  // v_BwdTrans(inarray,tmp);
843  //
844  // // weight = \grad \cdot tanvec
845  // for(int k = 0; k < GetCoordim(); ++k)
846  // {
847  // Vmath::Vcopy(nqtot, &(mkey.GetVariableCoefficient(0))[k*nqtot],
848  // 1, &tan[0], 1);
849  //
850  // // For Regular mesh ...
851  // if(gmatnumber==1)
852  // {
853  // // D_{/xi} and D_{/eta}
854  // v_PhysDeriv(0,tan,dtan0);
855  // v_PhysDeriv(1,tan,dtan1);
856  //
857  // // d v / d x_i = (d \xi / d x_i)*( d v / d \xi ) + (d \eta / d x_i)*( d v / d \eta )
858  // Vmath::Svtvp(nqtot,(mkey.GetVariableCoefficient(2*k+1))[0],&dtan0[0],1,&weight[0],1,&weight[0],1);
859  // Vmath::Svtvp(nqtot,(mkey.GetVariableCoefficient(2*k+2))[0],&dtan1[0],1,&weight[0],1,&weight[0],1);
860  // }
861  //
862  // // For Curved mesh ...
863  // else if(gmatnumber==nqtot)
864  // {
865  // // D_{x} and D_{y}
866  // v_PhysDeriv(k,tan,dtan0);
867  // Vmath::Vadd(nqtot,&dtan0[0],1,&weight[0],1,&weight[0],1);
868  // }
869  //
870  // else
871  // {
872  // ASSERTL1( ((gmatnumber=1) || (gmatnumber==nqtot) ), "Gmat is not in a right size");
873  // }
874  // }
875  //
876  // Vmath::Vmul(nqtot, &weight[0], 1, &tmp[0], 1, &tmp[0], 1);
877  // v_IProductWRTBase(tmp, outarray);
878  }
879 
881  Array<OneD,NekDouble> &outarray,
882  const StdMatrixKey &mkey,
883  bool addDiffusionTerm)
884  {
885 
886  int i;
887  int ndir = mkey.GetNVarCoeff(); // assume num.r consts corresponds to directions
888  ASSERTL0(ndir,"Must define at least one advection velocity");
889 
890  NekDouble lambda = mkey.GetConstFactor(eFactorLambda);
891  int totpts = GetTotPoints();
892  Array<OneD, NekDouble> tmp(3*totpts);
893  Array<OneD, NekDouble> tmp_deriv = tmp + totpts;
894  Array<OneD, NekDouble> tmp_adv = tmp_deriv + totpts;
895 
896 
897  ASSERTL1(ndir <= GetCoordim(),"Number of constants is larger than coordinate dimensions");
898 
899  v_BwdTrans(inarray,tmp);
900 
901  VarCoeffType varcoefftypes[] = {eVarCoeffVelX, eVarCoeffVelY};
902 
903  //calculate u dx + v dy + ..
904  Vmath::Zero(totpts,tmp_adv,1);
905  for(i = 0; i < ndir; ++i)
906  {
907  v_PhysDeriv(i,tmp,tmp_deriv);
908  Vmath::Vvtvp(totpts,mkey.GetVarCoeff(varcoefftypes[i]),1,tmp_deriv,1,tmp_adv,1,tmp_adv,1);
909  }
910 
911  if(lambda) // add -lambda*u
912  {
913  Vmath::Svtvp(totpts,-lambda,tmp,1,tmp_adv,1,tmp_adv,1);
914  }
915 
916 
917  if(addDiffusionTerm)
918  {
920  StdMatrixKey mkeylap(eLaplacian,DetShapeType(),*this);
921  LaplacianMatrixOp(inarray,lap,mkeylap);
922 
923  v_IProductWRTBase(tmp_adv, outarray);
924  // Lap v - u.grad v + lambda*u
925  // => (grad u, grad v) + u.grad v - lambda*u
926  Vmath::Vadd(m_ncoeffs,lap,1,outarray,1,outarray,1);
927  }
928  else
929  {
930  v_IProductWRTBase(tmp_adv, outarray);
931  }
932 
933  }
934 
935 
937  Array<OneD,NekDouble> &outarray,
938  const StdMatrixKey &mkey)
939  {
940  NekDouble lambda = mkey.GetConstFactor(eFactorLambda);
942  StdMatrixKey mkeymass(eMass,DetShapeType(),*this);
943  StdMatrixKey mkeylap(eLaplacian,DetShapeType(),*this);
944 
945  MassMatrixOp(inarray,tmp,mkeymass);
946  LaplacianMatrixOp(inarray,outarray,mkeylap);
947 
948  Blas::Daxpy(m_ncoeffs, lambda, tmp, 1, outarray, 1);
949  }
950 
952  Array<OneD, NekDouble> &outarray)
953  {
954  int nq = GetTotPoints();
955  StdMatrixKey bwdtransmatkey(eBwdTrans,DetShapeType(),*this);
956  DNekMatSharedPtr bwdtransmat = GetStdMatrix(bwdtransmatkey);
957 
958  Blas::Dgemv('N',nq,m_ncoeffs,1.0,bwdtransmat->GetPtr().get(),
959  nq, inarray.get(), 1, 0.0, outarray.get(), 1);
960  }
961 
962  // VIRTUAL INLINE FUNCTIONS FROM HEADER FILE
963  void StdExpansion::SetUpPhysNormals(const int edge)
964  {
965  v_SetUpPhysNormals(edge);
966  }
967 
969  const Array<OneD, const NekDouble> &physvals)
970  {
971  return v_StdPhysEvaluate(Lcoord,physvals);
972  }
973 
975  {
976  return m_elmt_id;
977  }
978 
980  {
981  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
982  return NullNekDouble1DArray;
983  }
984 
985 
987  {
988  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
989  }
990 
991  void StdExpansion::v_SetUpPhysNormals(const int edge)
992  {
993  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
994  }
995 
996  int StdExpansion::v_CalcNumberOfCoefficients(const std::vector<unsigned int> &nummodes, int &modes_offset)
997  {
998  NEKERROR(ErrorUtil::efatal, "This function is not defined for this class");
999  return 0;
1000  }
1001 
1003  const std::vector<unsigned int > &nummodes,
1004  const int nmode_offset,
1005  NekDouble *coeffs)
1006  {
1007  NEKERROR(ErrorUtil::efatal, "This function is not defined for this class");
1008  }
1009 
1011  {
1012  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
1013  }
1014 
1016  {
1017  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
1018  }
1019 
1021  const Array<OneD, const NekDouble> &Fy,
1022  const Array<OneD, const NekDouble> &Fz,
1023  Array< OneD, NekDouble> &outarray)
1024  {
1025  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
1026  }
1027 
1029  {
1030  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
1031  }
1032 
1033 
1035  {
1036  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
1038  }
1039 
1041  {
1042  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
1043  }
1044 
1046 
1047  {
1048  NEKERROR(ErrorUtil::efatal, "This function is only valid for three-dimensional LocalRegions");
1049  return eDir1FwdDir1_Dir2FwdDir2;
1050  }
1051 
1053  {
1054  NEKERROR(ErrorUtil::efatal, "This function is only valid for two-dimensional LocalRegions");
1055  return eForwards;
1056  }
1057 
1059  {
1060  NEKERROR(ErrorUtil::efatal, "This function is only valid for one-dimensional LocalRegions");
1061  return eFwd;
1062  }
1063 
1064 
1066  {
1067  NEKERROR(ErrorUtil::efatal, "This function is only valid for two-dimensional LocalRegions");
1068  return eForwards;
1069  }
1070 
1071 
1074  Array<OneD, NekDouble> &outarray)
1075  {
1076  NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1077  }
1078 
1080  Array<OneD, NekDouble> &coeffs,
1082  {
1083  NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1084  }
1085 
1086 
1088  const Array<OneD, const NekDouble> &physvals)
1089 
1090  {
1091  NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1092  return 0;
1093  }
1094 
1096  {
1097  NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
1098  }
1099 
1101  {
1102  ASSERTL0(false, "This function is needs defining for this shape");
1103  return 0;
1104  }
1105 
1107  {
1108  ASSERTL0(false, "This function is needs defining for this shape");
1109  return 0;
1110  }
1111 
1113  {
1114  ASSERTL0(false, "This function is needs defining for this shape");
1115  return 0;
1116  }
1117 
1119  {
1120  ASSERTL0(false, "This function is needs defining for this shape");
1121  return 0;
1122  }
1123 
1124  int StdExpansion::v_GetEdgeNcoeffs(const int i) const
1125  {
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  ASSERTL0(false, "This function is not valid or not defined");
1139  return 0;
1140  }
1141 
1143  {
1144  ASSERTL0(false, "This function is not valid or not defined");
1145  return 0;
1146  }
1147 
1149  {
1150  ASSERTL0(false, "This function is not valid or not defined");
1152  }
1153 
1154  const LibUtilities::BasisKey StdExpansion::v_DetFaceBasisKey(const int i, const int k) const
1155  {
1156  ASSERTL0(false, "This function is not valid or not defined");
1158  }
1159 
1160  int StdExpansion::v_GetFaceNumPoints(const int i) const
1161  {
1162  ASSERTL0(false, "This function is not valid or not defined");
1163  return 0;
1164  }
1165 
1166  int StdExpansion::v_GetFaceNcoeffs(const int i) const
1167  {
1168  ASSERTL0(false, "This function is not valid or not defined");
1169  return 0;
1170  }
1171 
1172  int StdExpansion::v_GetFaceIntNcoeffs(const int i) const
1173  {
1174  ASSERTL0(false, "This function is not valid or not defined");
1175  return 0;
1176  }
1177 
1179  {
1180  ASSERTL0(false, "This function is not valid or not defined");
1181  return 0;
1182  }
1183 
1184  int StdExpansion::v_GetTraceNcoeffs(const int i) const
1185  {
1186  ASSERTL0(false, "This function is not valid or not defined");
1187  return 0;
1188  }
1189 
1190 
1192  {
1193  ASSERTL0(false, "This function is not valid or not defined");
1195  }
1196 
1198  {
1199  ASSERTL0(false, "This function is not valid or not defined");
1200 
1202  }
1203 
1205  {
1206  ASSERTL0(false, "This function is not valid or not defined");
1207 
1209  }
1210 
1212  {
1213  ASSERTL0(false, "This expansion does not have a shape type defined");
1215  }
1216 
1217  boost::shared_ptr<StdExpansion>
1219  {
1220  ASSERTL0(false,"This method is not defined for this expansion");
1221  StdExpansionSharedPtr returnval;
1222  return returnval;
1223  }
1224 
1226  {
1227  ASSERTL0(false, "This function is not valid or not defined");
1228  return 0;
1229  }
1230 
1232  {
1233  ASSERTL0(false,"This function has not been defined for this expansion");
1234  return false;
1235  }
1236 
1237 
1239  {
1240  return false;
1241  }
1242 
1244  const Array<OneD, const NekDouble>& inarray,
1245  Array<OneD, NekDouble> &outarray)
1246  {
1247  NEKERROR(ErrorUtil::efatal, "This method has not been defined");
1248  }
1249 
1250  /**
1251  *
1252  */
1254  Array<OneD, NekDouble> &outarray)
1255  {
1256  NEKERROR(ErrorUtil::efatal, "This method has not been defined");
1257  }
1258 
1259 
1260  /**
1261  * @brief Integrates the specified function over the domain.
1262  * @see StdRegions#StdExpansion#Integral.
1263  */
1265  {
1266  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1267  "local expansions");
1268  return 0;
1269  }
1270 
1271 
1272  void StdExpansion::v_AddRobinMassMatrix(const int edgeid, const Array<OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
1273  {
1274  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1275  "specific element types");
1276  }
1277 
1279  {
1280  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1281  "specific element types");
1282  }
1283 
1284  /**
1285  * @brief Calculate the derivative of the physical points
1286  * @see StdRegions#StdExpansion#PhysDeriv
1287  */
1289  Array<OneD, NekDouble> &out_d1,
1290  Array<OneD, NekDouble> &out_d2,
1291  Array<OneD, NekDouble> &out_d3)
1292  {
1293  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1294  "local expansions");
1295  }
1296 
1298  Array<OneD, NekDouble> &out_ds)
1299  {
1300  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1301  "local expansions");
1302  }
1304  Array<OneD, NekDouble>& out_dn)
1305  {
1306  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1307  "local expansions");
1308  }
1309 
1310  /**
1311  * @brief Calculate the derivative of the physical points in a
1312  * given direction
1313  * @see StdRegions#StdExpansion#PhysDeriv
1314  */
1315  void StdExpansion::v_PhysDeriv(const int dir,
1316  const Array<OneD, const NekDouble>& inarray,
1317  Array<OneD, NekDouble> &out_d0)
1318 
1319  {
1320  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1321  "specific element types");
1322  }
1323 
1324  /**
1325  * @brief Physical derivative along a direction vector.
1326  * @see StdRegions#StdExpansion#PhysDirectionalDeriv
1327  */
1329  const Array<OneD, const NekDouble>& direction,
1330  Array<OneD, NekDouble> &outarray)
1331  {
1332  NEKERROR(ErrorUtil::efatal, "This function is only valid for "
1333  "specific element types");
1334  }
1335 
1337  Array<OneD, NekDouble> &out_d1,
1338  Array<OneD, NekDouble> &out_d2,
1339  Array<OneD, NekDouble> &out_d3)
1340  {
1341  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1342  }
1343 
1344  void StdExpansion::v_StdPhysDeriv (const int dir,
1345  const Array<OneD, const NekDouble>& inarray,
1346  Array<OneD, NekDouble> &outarray)
1347  {
1348  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1349  }
1350 
1351 
1353  {
1354  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1355  return 0;
1356  }
1357 
1358 
1360  {
1361  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1362  return 0;
1363  }
1364 
1365 
1366  void StdExpansion::v_FillMode(const int mode, Array<OneD, NekDouble> &outarray)
1367  {
1368  NEKERROR(ErrorUtil::efatal, "This function has not "
1369  "been defined for this shape");
1370  }
1371 
1373  {
1374  NEKERROR(ErrorUtil::efatal, "This function has not "
1375  "been defined for this element");
1376  DNekMatSharedPtr returnval;
1377  return returnval;
1378  }
1379 
1381  {
1382  NEKERROR(ErrorUtil::efatal, "This function has not "
1383  "been defined for this element");
1384  DNekMatSharedPtr returnval;
1385  return returnval;
1386  }
1387 
1389  Array<OneD, NekDouble> &coords_1,
1390  Array<OneD, NekDouble> &coords_2)
1391  {
1392  NEKERROR(ErrorUtil::efatal, "Write coordinate definition method");
1393  }
1394 
1396  Array<OneD, NekDouble> &coord)
1397  {
1398  NEKERROR(ErrorUtil::efatal, "Write coordinate definition method");
1399  }
1400 
1402  {
1403  NEKERROR(ErrorUtil::efatal, "Write method");
1404  return -1;
1405  }
1406 
1408  {
1409  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1410  }
1411 
1413  {
1414  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1415  }
1416 
1417  int StdExpansion::v_GetVertexMap(const int localVertexId,
1418  bool useCoeffPacking)
1419  {
1420  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1421  return 0;
1422  }
1423 
1424  void StdExpansion::v_GetEdgeInteriorMap(const int eid, const Orientation edgeOrient,
1425  Array<OneD, unsigned int> &maparray,
1426  Array<OneD, int> &signarray)
1427  {
1428  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1429  }
1430 
1431  void StdExpansion::v_GetFaceInteriorMap(const int fid, const Orientation faceOrient,
1432  Array<OneD, unsigned int> &maparray,
1433  Array<OneD, int> &signarray)
1434  {
1435  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1436  }
1437 
1439  const int eid,
1440  const Orientation edgeOrient,
1441  Array<OneD, unsigned int>& maparray,
1442  Array<OneD, int>& signarray,
1443  int P)
1444  {
1445  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape");
1446  }
1447 
1448  void StdExpansion::v_GetFaceToElementMap(const int fid, const Orientation faceOrient,
1449  Array<OneD, unsigned int> &maparray,
1450  Array<OneD, int> &signarray,
1451  int nummodesA, int nummodesB)
1452  {
1453  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1454  }
1455 
1457  {
1458  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape or library" );
1459  }
1460 
1461  void StdExpansion::v_GetEdgePhysVals(const int edge, const boost::shared_ptr<StdExpansion> &EdgeExp, const Array<OneD, const NekDouble> &inarray, Array<OneD,NekDouble> &outarray)
1462  {
1463  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape or library" );
1464  }
1465 
1466  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)
1467  {
1468  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape or library" );
1469  }
1470 
1471  void StdExpansion::v_GetVertexPhysVals(const int vertex, const Array<OneD, const NekDouble> &inarray, NekDouble &outarray)
1472  {
1473  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape or library" );
1474  }
1475 
1477  {
1478  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape or library" );
1479  }
1480 
1482  const int edge,
1483  Array<OneD, NekDouble> &outarray)
1484  {
1486  "Method does not exist for this shape or library");
1487  }
1488 
1489  void StdExpansion::v_GetFacePhysVals( const int face,
1490  const boost::shared_ptr<StdExpansion> &FaceExp,
1491  const Array<OneD, const NekDouble> &inarray,
1492  Array<OneD, NekDouble> &outarray,
1493  StdRegions::Orientation orient)
1494  {
1495  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape or library" );
1496  }
1497 
1499  const int edge,
1500  Array<OneD, int> &outarray)
1501  {
1503  "Method does not exist for this shape or library" );
1504  }
1505 
1506  void StdExpansion::v_GetFacePhysMap(const int face,
1507  Array<OneD, int> &outarray)
1508  {
1509  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape or library" );
1510  }
1511 
1513  const Array<OneD, const NekDouble> &inarray,
1514  Array<OneD, NekDouble> &outarray)
1515  {
1516  v_MultiplyByStdQuadratureMetric(inarray,outarray);
1517  }
1518 
1520  const Array<OneD, const NekDouble> &inarray,
1521  Array<OneD, NekDouble> &outarray)
1522  {
1523  NEKERROR(ErrorUtil::efatal, "Method does not exist for this shape or library");
1524  }
1525 
1526  const boost::shared_ptr<SpatialDomains::GeomFactors>& StdExpansion::v_GetMetricInfo() const
1527  {
1528  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
1530 
1531  }
1532 
1534  Array<OneD, NekDouble> &outarray)
1535  {
1536  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1537  }
1538 
1540  Array<OneD, NekDouble> &outarray,
1541  bool multiplybyweights)
1542  {
1543  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1544  }
1545 
1547  const Array<OneD, const NekDouble>& inarray,
1548  Array<OneD, NekDouble> &outarray)
1549  {
1550  NEKERROR(ErrorUtil::efatal,"Method does not exist for this shape" );
1551  }
1552 
1554  Array<OneD,NekDouble> &outarray,
1555  const StdMatrixKey &mkey)
1556  {
1557  // If this function is not reimplemented on shape level, the function
1558  // below will be called
1559  MassMatrixOp_MatFree(inarray,outarray,mkey);
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  LaplacianMatrixOp_MatFree(inarray,outarray,mkey);
1569  }
1570 
1571 
1573  const StdMatrixKey &mkey)
1574  {
1575  ASSERTL0(false, "This function is not defined in StdExpansion.");
1576  }
1577 
1579  const Array<OneD, const NekDouble> &inarray,
1580  Array<OneD, NekDouble> &outarray)
1581  {
1582  ASSERTL0(false, "This function is not defined in StdExpansion.");
1583  }
1584 
1585  void StdExpansion::v_LaplacianMatrixOp(const int k1, const int k2,
1586  const Array<OneD, const NekDouble> &inarray,
1587  Array<OneD,NekDouble> &outarray,
1588  const StdMatrixKey &mkey)
1589  {
1590  // If this function is not reimplemented on shape level, the function
1591  // below will be called
1592  LaplacianMatrixOp_MatFree(k1,k2,inarray,outarray,mkey);
1593  }
1594 
1596  const Array<OneD, const NekDouble> &inarray,
1597  Array<OneD,NekDouble> &outarray,
1598  const StdMatrixKey &mkey)
1599  {
1600  // If this function is not reimplemented on shape level, the function
1601  // below will be called
1602  WeakDerivMatrixOp_MatFree(i,inarray,outarray,mkey);
1603 
1604  }
1605 
1607  Array<OneD,NekDouble> &outarray,
1608  const StdMatrixKey &mkey)
1609  {
1610  // If this function is not reimplemented on shape level, the function
1611  // below will be called
1612  WeakDirectionalDerivMatrixOp_MatFree(inarray,outarray,mkey);
1613 
1614  }
1615 
1617  Array<OneD,NekDouble> &outarray,
1618  const StdMatrixKey &mkey)
1619  {
1620  // If this function is not reimplemented on shape level, the function
1621  // below will be called
1622  MassLevelCurvatureMatrixOp_MatFree(inarray,outarray,mkey);
1623  }
1624 
1626  const NekDouble> &inarray,
1627  Array<OneD,NekDouble> &outarray,
1628  const StdMatrixKey &mkey, bool addDiffusionTerm)
1629  {
1630  // If this function is not reimplemented on shape level, the function
1631  // below will be called
1632  LinearAdvectionDiffusionReactionMatrixOp_MatFree(inarray,outarray,mkey,addDiffusionTerm);
1633 
1634  }
1635 
1637  Array<OneD,NekDouble> &outarray,
1638  const StdMatrixKey &mkey)
1639  {
1640  // If this function is not reimplemented on shape level, the function
1641  // below will be called
1642  HelmholtzMatrixOp_MatFree(inarray,outarray,mkey);
1643  }
1644 
1646  Array<OneD,NekDouble> &outarray,
1647  const StdMatrixKey &mkey)
1648  {
1649  // If this function is not reimplemented on shape level, the function
1650  // below will be called
1651  LaplacianMatrixOp_MatFree_GenericImpl(inarray,outarray,mkey);
1652  }
1653 
1655  const Array<OneD, const NekDouble> &inarray,
1656  Array<OneD, NekDouble> &outarray,
1658  {
1659  ASSERTL0(false, "Not implemented.");
1660  }
1661 
1663  Array<OneD,NekDouble> &outarray,
1664  const StdMatrixKey &mkey)
1665  {
1666  // If this function is not reimplemented on shape level, the function
1667  // below will be called
1668  HelmholtzMatrixOp_MatFree_GenericImpl(inarray,outarray,mkey);
1669  }
1670 
1671  const NormalVector & StdExpansion::v_GetEdgeNormal(const int edge) const
1672  {
1673  ASSERTL0(false, "Cannot get edge normals for this expansion.");
1674  static NormalVector result;
1675  return result;
1676  }
1677 
1679  {
1680  ASSERTL0(false, "Cannot compute edge normal for this expansion.");
1681  }
1682 
1684  {
1685  ASSERTL0(false, "Not implemented.");
1686  }
1687 
1689  {
1690  ASSERTL0(false, "Not implemented.");
1691  return false;
1692  }
1693 
1695  {
1696  ASSERTL0(false, "Cannot compute face normal for this expansion.");
1697  }
1698 
1700  {
1701  ASSERTL0(false, "Not implemented.");
1702  }
1703 
1705  {
1706  ASSERTL0(false, "Not implemented.");
1707  return false;
1708  }
1709 
1711  {
1712  ASSERTL0(false, "Cannot compute vertex normal for this expansion.");
1713  }
1714 
1715  const NormalVector & StdExpansion::v_GetFaceNormal(const int face) const
1716  {
1717  ASSERTL0(false, "Cannot get face normals for this expansion.");
1718  static NormalVector result;
1719  return result;
1720  }
1721 
1722  const NormalVector & StdExpansion::v_GetVertexNormal(const int vertex) const
1723  {
1724  ASSERTL0(false, "Cannot get vertex normals for this expansion.");
1725  static NormalVector result;
1726  return result;
1727  }
1728 
1730  {
1731  ASSERTL0(false, "Cannot get face normals for this expansion.");
1732  static NormalVector result;
1733  return result;
1734  }
1735 
1738  {
1739  ASSERTL0(false, "Not implemented.");
1740  Array<OneD, unsigned int> noinversemap(1);
1741  return noinversemap;
1742  }
1743 
1746  StdRegions::Orientation faceOrient)
1747  {
1748  ASSERTL0(false, "Not implemented.");
1749  Array<OneD, unsigned int> noinversemap(1);
1750  return noinversemap;
1751  }
1752 
1755  const DNekScalMatSharedPtr & m_transformationmatrix)
1756  {
1757  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
1758  return NullDNekMatSharedPtr;
1759  }
1760 
1762  const Array<OneD, const NekDouble> &inarray,
1763  Array<OneD, NekDouble> &outarray,
1764  int npset)
1765  {
1767  DNekMatSharedPtr intmat;
1768 
1769  int nqtot = GetTotPoints();
1770  int np = 0;
1771  if(npset == -1) // use values from basis num points()
1772  {
1773  int nqbase;
1774  for(int i = 0; i < m_base.num_elements(); ++i)
1775  {
1776  nqbase = m_base[i]->GetNumPoints();
1777  np = std::max(np,nqbase);
1778  }
1779 
1780  StdMatrixKey Ikey(ePhysInterpToEquiSpaced, shape, *this);
1781  intmat = GetStdMatrix(Ikey);
1782  }
1783  else
1784  {
1785  np = npset;
1786 
1787  ConstFactorMap cmap;
1788  cmap[eFactorConst] = np;
1789  StdMatrixKey Ikey(ePhysInterpToEquiSpaced, shape, *this, cmap);
1790  intmat = GetStdMatrix(Ikey);
1791 
1792  }
1793 
1794  NekVector<NekDouble> in (nqtot,inarray,eWrapper);
1796  out = (*intmat) * in;
1797  }
1798 
1800  Array<OneD, int> &conn,
1801  bool standard)
1802  {
1803  ASSERTL0(false, "Not implemented.");
1804  }
1805 
1807  const Array<OneD, const NekDouble> &inarray,
1808  Array<OneD, NekDouble> &outarray)
1809  {
1811 
1812  // inarray has to be consistent with NumModes definition
1813  // There is also a check in GetStdMatrix to see if all
1814  // modes are of the same size
1815  ConstFactorMap cmap;
1816 
1817  cmap[eFactorConst] = m_base[0]->GetNumModes();
1818  StdMatrixKey Ikey(eEquiSpacedToCoeffs, shape, *this,cmap);
1819  DNekMatSharedPtr intmat = GetStdMatrix(Ikey);
1820 
1821  NekVector<NekDouble> in (m_ncoeffs, inarray, eWrapper);
1822  NekVector<NekDouble> out(m_ncoeffs, outarray,eWrapper);
1823  out = (*intmat) * in;
1824  }
1825 
1826  }//end namespace
1827 }//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:470
NekDouble GetConstFactor(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:122
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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 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:592
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:971
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:992
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.
virtual void v_GetEdgeToElementMap(const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int P=-1)
LibUtilities::ShapeType GetShapeType() const
Definition: StdMatrixKey.h:87
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:251
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)
virtual void v_NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
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:700
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
void IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
int GetNumberOfCoefficients(ShapeType shape, std::vector< unsigned int > &modes, int offset)
Definition: ShapeType.hpp:312
virtual void v_GetEdgePhysMap(const int edge, Array< OneD, int > &outarray)
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)
void MassLevelCurvatureMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
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:821
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:811
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
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...
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:978
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 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)
#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)
virtual bool v_FaceNormalNegated(const int face)
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)
virtual int v_GetNfaces() const
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:359
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...
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:1047
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:816
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:227
virtual int v_GetTraceNcoeffs(const int i) const
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:252
void LinearAdvectionDiffusionReactionMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
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)