Nektar++
StdQuadExp.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File StdQuadExp.cpp
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: Quadrilateral routines built upon StdExpansion2D
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <boost/core/ignore_unused.hpp>
36 
38 #include <StdRegions/StdQuadExp.h>
39 #include <StdRegions/StdSegExp.h>
41 
42 using namespace std;
43 
44 namespace Nektar
45 {
46  namespace StdRegions
47  {
48 
49 
50  StdQuadExp::StdQuadExp()
51  {
52  }
53 
54  /** \brief Constructor using BasisKey class for quadrature
55  * points and order definition
56  */
57  StdQuadExp::StdQuadExp(const LibUtilities::BasisKey &Ba,
58  const LibUtilities::BasisKey &Bb):
59  StdExpansion (Ba.GetNumModes()*Bb.GetNumModes(),2,Ba,Bb),
60  StdExpansion2D(Ba.GetNumModes()*Bb.GetNumModes(),Ba,Bb)
61  {
62  }
63 
64  /** \brief Copy Constructor */
66  StdExpansion(T),
68  {
69  }
70 
71  /** \brief Destructor */
73  {
74  }
75 
76  /////////////////////////
77  // Integration Methods //
78  /////////////////////////
79 
81  const Array<OneD, const NekDouble>& inarray)
82  {
83  Array<OneD, const NekDouble> w0 = m_base[0]->GetW();
84  Array<OneD, const NekDouble> w1 = m_base[1]->GetW();
85 
86  return StdExpansion2D::Integral(inarray,w0,w1);
87  }
88 
89  /////////////////////////////
90  // Differentiation Methods //
91  /////////////////////////////
92 
93  /** \brief Calculate the derivative of the physical points
94  *
95  * For quadrilateral region can use the Tensor_Deriv function
96  * defined under StdExpansion.
97  */
98 
100  Array<OneD, NekDouble> &out_d0,
101  Array<OneD, NekDouble> &out_d1,
102  Array<OneD, NekDouble> &out_d2)
103  {
104  boost::ignore_unused(out_d2);
105  PhysTensorDeriv(inarray, out_d0, out_d1);
106  }
107 
108  void StdQuadExp::v_PhysDeriv(const int dir,
109  const Array<OneD, const NekDouble>& inarray,
110  Array<OneD, NekDouble> &outarray)
111  {
112  switch(dir)
113  {
114  case 0:
115  {
116  PhysTensorDeriv(inarray, outarray, NullNekDouble1DArray);
117  }
118  break;
119  case 1:
120  {
121  PhysTensorDeriv(inarray, NullNekDouble1DArray, outarray);
122  }
123  break;
124  default:
125  {
126  ASSERTL1(false,"input dir is out of range");
127  }
128  break;
129  }
130  }
131 
133  Array<OneD, NekDouble> &out_d0,
134  Array<OneD, NekDouble> &out_d1,
135  Array<OneD, NekDouble> &out_d2)
136  {
137  boost::ignore_unused(out_d2);
138  //PhysTensorDeriv(inarray, out_d0, out_d1);
139  StdQuadExp::v_PhysDeriv(inarray, out_d0, out_d1);
140  }
141 
142  void StdQuadExp::v_StdPhysDeriv(const int dir,
143  const Array<OneD, const NekDouble>& inarray,
144  Array<OneD, NekDouble> &outarray)
145  {
146  //PhysTensorDeriv(inarray, outarray);
147  StdQuadExp::v_PhysDeriv(dir,inarray,outarray);
148 
149  }
150 
151 
152  ////////////////
153  // Transforms //
154  ////////////////
155 
157  Array<OneD, NekDouble> &outarray)
158  {
159  if(m_base[0]->Collocation() && m_base[1]->Collocation())
160  {
162  inarray, 1, outarray, 1);
163  }
164  else
165  {
166  StdQuadExp::v_BwdTrans_SumFac(inarray,outarray);
167  }
168  }
169 
171  const Array<OneD, const NekDouble>& inarray,
172  Array<OneD, NekDouble> &outarray)
173  {
175  m_base[1]->GetNumModes());
176 
177  BwdTrans_SumFacKernel(m_base[0]->GetBdata(),
178  m_base[1]->GetBdata(),
179  inarray,outarray,wsp,true,true);
180  }
181 
182  // The arguments doCheckCollDir0 and doCheckCollDir1 allow you to specify whether
183  // to check if the basis has collocation properties (i.e. for the classical spectral
184  // element basis, In this case the 1D 'B' matrix is equal to the identity matrix
185  // which can be exploited to speed up the calculations).
186  // However, as this routine also allows to pass the matrix 'DB' (derivative of the basis),
187  // the collocation property cannot always be used. Therefor follow this rule:
188  // if base0 == m_base[0]->GetBdata() --> set doCheckCollDir0 == true;
189  // base1 == m_base[1]->GetBdata() --> set doCheckCollDir1 == true;
190  // base0 == m_base[0]->GetDbdata() --> set doCheckCollDir0 == false;
191  // base1 == m_base[1]->GetDbdata() --> set doCheckCollDir1 == false;
193  const Array<OneD, const NekDouble>& base0,
194  const Array<OneD, const NekDouble>& base1,
195  const Array<OneD, const NekDouble>& inarray,
196  Array<OneD, NekDouble> &outarray,
198  bool doCheckCollDir0,
199  bool doCheckCollDir1)
200  {
201  int nquad0 = m_base[0]->GetNumPoints();
202  int nquad1 = m_base[1]->GetNumPoints();
203  int nmodes0 = m_base[0]->GetNumModes();
204  int nmodes1 = m_base[1]->GetNumModes();
205 
206  bool colldir0 = doCheckCollDir0?(m_base[0]->Collocation()):false;
207  bool colldir1 = doCheckCollDir1?(m_base[1]->Collocation()):false;
208 
209  if(colldir0 && colldir1)
210  {
211  Vmath::Vcopy(m_ncoeffs,inarray.get(),1,outarray.get(),1);
212  }
213  else if(colldir0)
214  {
215  Blas::Dgemm('N','T', nquad0, nquad1,nmodes1, 1.0, &inarray[0], nquad0,
216  base1.get(), nquad1, 0.0, &outarray[0], nquad0);
217  }
218  else if(colldir1)
219  {
220  Blas::Dgemm('N','N', nquad0,nmodes1,nmodes0,1.0, base0.get(),
221  nquad0, &inarray[0], nmodes0,0.0,&outarray[0], nquad0);
222  }
223  else
224  {
225  ASSERTL1(wsp.num_elements()>=nquad0*nmodes1,"Workspace size is not sufficient");
226 
227  // Those two calls correpsond to the operation
228  // out = B0*in*Transpose(B1);
229  Blas::Dgemm('N','N', nquad0,nmodes1,nmodes0,1.0, base0.get(),
230  nquad0, &inarray[0], nmodes0,0.0,&wsp[0], nquad0);
231  Blas::Dgemm('N','T', nquad0, nquad1,nmodes1, 1.0, &wsp[0], nquad0,
232  base1.get(), nquad1, 0.0, &outarray[0], nquad0);
233  }
234  }
235 
236 
238  Array<OneD, NekDouble> &outarray)
239  {
240  if((m_base[0]->Collocation())&&(m_base[1]->Collocation()))
241  {
242  Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
243  }
244  else
245  {
246  StdQuadExp::v_IProductWRTBase(inarray,outarray);
247 
248  // get Mass matrix inverse
249  StdMatrixKey masskey(eInvMass,DetShapeType(),*this);
250  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
251 
252  // copy inarray in case inarray == outarray
253  NekVector<NekDouble> in(m_ncoeffs,outarray,eCopy);
255 
256  out = (*matsys)*in;
257  }
258  }
259 
261  const Array<OneD, const NekDouble>& inarray,
262  Array<OneD, NekDouble> &outarray)
263  {
264  if((m_base[0]->Collocation())&&(m_base[1]->Collocation()))
265  {
266  Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
267  }
268  else
269  {
270  int i,j;
271  int npoints[2] = {m_base[0]->GetNumPoints(),
272  m_base[1]->GetNumPoints()};
273  int nmodes[2] = {m_base[0]->GetNumModes(),
274  m_base[1]->GetNumModes()};
275 
276  fill(outarray.get(), outarray.get()+m_ncoeffs, 0.0 );
277 
278  Array<OneD, NekDouble> physEdge[4];
279  Array<OneD, NekDouble> coeffEdge[4];
280  for(i = 0; i < 4; i++)
281  {
282  physEdge[i] = Array<OneD, NekDouble>(npoints[i%2]);
283  coeffEdge[i] = Array<OneD, NekDouble>(nmodes[i%2]);
284  }
285 
286  for(i = 0; i < npoints[0]; i++)
287  {
288  physEdge[0][i] = inarray[i];
289  physEdge[2][i] = inarray[npoints[0]*npoints[1]-1-i];
290  }
291 
292  for(i = 0; i < npoints[1]; i++)
293  {
294  physEdge[1][i] = inarray[npoints[0]-1+i*npoints[0]];
295  physEdge[3][i] = inarray[(npoints[1]-1)*npoints[0]-i*npoints[0]];
296  }
297 
300 
301  Array<OneD, unsigned int> mapArray;
302  Array<OneD, int> signArray;
303  NekDouble sign;
304 
305  for(i = 0; i < 4; i++)
306  {
307  segexp[i%2]->FwdTrans_BndConstrained(physEdge[i],coeffEdge[i]);
308 
309  GetEdgeToElementMap(i,eForwards,mapArray,signArray);
310  for(j=0; j < nmodes[i%2]; j++)
311  {
312  sign = (NekDouble) signArray[j];
313  outarray[ mapArray[j] ] = sign * coeffEdge[i][j];
314  }
315  }
316 
319 
320  StdMatrixKey masskey(eMass,DetShapeType(),*this);
321  MassMatrixOp(outarray,tmp0,masskey);
322  IProductWRTBase(inarray,tmp1);
323 
324  Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
325 
326  // get Mass matrix inverse (only of interior DOF)
327  // use block (1,1) of the static condensed system
328  // note: this block alreay contains the inverse matrix
329  DNekMatSharedPtr matsys = (m_stdStaticCondMatrixManager[masskey])->GetBlock(1,1);
330 
331  int nBoundaryDofs = NumBndryCoeffs();
332  int nInteriorDofs = m_ncoeffs - nBoundaryDofs;
333 
334 
335  Array<OneD, NekDouble> rhs(nInteriorDofs);
336  Array<OneD, NekDouble> result(nInteriorDofs);
337 
338  GetInteriorMap(mapArray);
339 
340  for(i = 0; i < nInteriorDofs; i++)
341  {
342  rhs[i] = tmp1[ mapArray[i] ];
343  }
344 
345  Blas::Dgemv('N',nInteriorDofs,nInteriorDofs,1.0, &(matsys->GetPtr())[0],
346  nInteriorDofs,rhs.get(),1,0.0,result.get(),1);
347 
348  for(i = 0; i < nInteriorDofs; i++)
349  {
350  outarray[ mapArray[i] ] = result[i];
351  }
352  }
353 
354  }
355 
356  /////////////////////////////
357  // Inner Product Functions //
358  /////////////////////////////
359 
360  /** \brief Calculate the inner product of inarray with respect to
361  * the basis B=base0*base1 and put into outarray
362  *
363  * \f$
364  * \begin{array}{rcl}
365  * I_{pq} = (\phi_q \phi_q, u) & = & \sum_{i=0}^{nq_0}
366  * \sum_{j=0}^{nq_1}
367  * \phi_p(\xi_{0,i}) \phi_q(\xi_{1,j}) w^0_i w^1_j u(\xi_{0,i}
368  * \xi_{1,j}) \\
369  * & = & \sum_{i=0}^{nq_0} \phi_p(\xi_{0,i})
370  * \sum_{j=0}^{nq_1} \phi_q(\xi_{1,j}) \tilde{u}_{i,j}
371  * \end{array}
372  * \f$
373  *
374  * where
375  *
376  * \f$ \tilde{u}_{i,j} = w^0_i w^1_j u(\xi_{0,i},\xi_{1,j}) \f$
377  *
378  * which can be implemented as
379  *
380  * \f$ f_{qi} = \sum_{j=0}^{nq_1} \phi_q(\xi_{1,j})
381  * \tilde{u}_{i,j} = {\bf B_1 U} \f$
382  * \f$ I_{pq} = \sum_{i=0}^{nq_0} \phi_p(\xi_{0,i}) f_{qi} =
383  * {\bf B_0 F} \f$
384  */
386  const Array<OneD, const NekDouble>& inarray,
387  Array<OneD, NekDouble> &outarray)
388  {
389  if(m_base[0]->Collocation() && m_base[1]->Collocation())
390  {
391  MultiplyByQuadratureMetric(inarray,outarray);
392  }
393  else
394  {
395  StdQuadExp::v_IProductWRTBase_SumFac(inarray,outarray);
396  }
397  }
398 
400  const Array<OneD, const NekDouble>& inarray,
401  Array<OneD, NekDouble> &outarray,
402  bool multiplybyweights)
403  {
404  int nquad0 = m_base[0]->GetNumPoints();
405  int nquad1 = m_base[1]->GetNumPoints();
406  int order0 = m_base[0]->GetNumModes();
407 
408  if(multiplybyweights)
409  {
410  Array<OneD,NekDouble> tmp(nquad0*nquad1+nquad1*order0);
411  Array<OneD,NekDouble> wsp(tmp+nquad0*nquad1);
412 
413  // multiply by integration constants
414  MultiplyByQuadratureMetric(inarray,tmp);
415  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
416  m_base[1]->GetBdata(),
417  tmp,outarray,wsp,true,true);
418  }
419  else
420  {
421  Array<OneD,NekDouble> wsp(nquad1*order0);
422  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
423  m_base[1]->GetBdata(),
424  inarray,outarray,wsp,true,true);
425  }
426  }
427 
429  const Array<OneD, const NekDouble>& inarray,
430  Array<OneD, NekDouble> &outarray)
431  {
432  int nq = GetTotPoints();
433  StdMatrixKey iprodmatkey(eIProductWRTBase,DetShapeType(),*this);
434  DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
435 
436  Blas::Dgemv('N',m_ncoeffs,nq,1.0,iprodmat->GetPtr().get(),
437  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
438  }
439 
441  const Array<OneD, const NekDouble>& inarray,
442  Array<OneD, NekDouble> & outarray)
443  {
444  StdQuadExp::IProductWRTDerivBase_SumFac(dir,inarray,outarray);
445  }
446 
448  const Array<OneD, const NekDouble>& inarray,
449  Array<OneD, NekDouble> &outarray)
450  {
451  ASSERTL0((dir==0)||(dir==1),"input dir is out of range");
452 
453  int nquad0 = m_base[0]->GetNumPoints();
454  int nquad1 = m_base[1]->GetNumPoints();
455  int nqtot = nquad0*nquad1;
456  int order0 = m_base[0]->GetNumModes();
457 
458  Array<OneD,NekDouble> tmp(nqtot+nquad1*order0);
459  Array<OneD,NekDouble> wsp(tmp+nqtot);
460 
461  // multiply by integration constants
462  MultiplyByQuadratureMetric(inarray,tmp);
463 
464  if(dir) // dir == 1
465  {
466  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
467  m_base[1]->GetDbdata(),
468  tmp,outarray,wsp,true,false);
469  }
470  else // dir == 0
471  {
472  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
473  m_base[1]->GetBdata(),
474  tmp,outarray,wsp,false,true);
475  }
476  }
477 
479  const Array<OneD, const NekDouble>& inarray,
480  Array<OneD, NekDouble> &outarray)
481  {
482  ASSERTL0((dir==0)||(dir==1),"input dir is out of range");
483 
484  int nq = GetTotPoints();
485  MatrixType mtype;
486 
487  if(dir) // dir == 1
488  {
489  mtype = eIProductWRTDerivBase1;
490  }
491  else // dir == 0
492  {
493  mtype = eIProductWRTDerivBase0;
494  }
495 
496  StdMatrixKey iprodmatkey(mtype,DetShapeType(),*this);
497  DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
498 
499  Blas::Dgemv('N',m_ncoeffs,nq,1.0,iprodmat->GetPtr().get(),
500  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
501  }
502 
503  // the arguments doCheckCollDir0 and doCheckCollDir1 allow you to specify whether
504  // to check if the basis has collocation properties (i.e. for the classical spectral
505  // element basis, In this case the 1D 'B' matrix is equal to the identity matrix
506  // which can be exploited to speed up the calculations).
507  // However, as this routine also allows to pass the matrix 'DB' (derivative of the basis),
508  // the collocation property cannot always be used. Therefor follow this rule:
509  // if base0 == m_base[0]->GetBdata() --> set doCheckCollDir0 == true;
510  // base1 == m_base[1]->GetBdata() --> set doCheckCollDir1 == true;
511  // base0 == m_base[0]->GetDbdata() --> set doCheckCollDir0 == false;
512  // base1 == m_base[1]->GetDbdata() --> set doCheckCollDir1 == false;
514  const Array<OneD, const NekDouble>& base0,
515  const Array<OneD, const NekDouble>& base1,
516  const Array<OneD, const NekDouble>& inarray,
517  Array<OneD, NekDouble> &outarray,
519  bool doCheckCollDir0,
520  bool doCheckCollDir1)
521  {
522  int nquad0 = m_base[0]->GetNumPoints();
523  int nquad1 = m_base[1]->GetNumPoints();
524  int nmodes0 = m_base[0]->GetNumModes();
525  int nmodes1 = m_base[1]->GetNumModes();
526 
527  bool colldir0 = doCheckCollDir0?(m_base[0]->Collocation()):false;
528  bool colldir1 = doCheckCollDir1?(m_base[1]->Collocation()):false;
529 
530  if(colldir0 && colldir1)
531  {
532  Vmath::Vcopy(m_ncoeffs,inarray.get(),1,outarray.get(),1);
533  }
534  else if(colldir0)
535  {
536  Blas::Dgemm('N', 'N',nmodes0,nmodes1, nquad1,1.0, inarray.get(),
537  nmodes0, base1.get(), nquad1, 0.0,outarray.get(),nmodes0);
538  }
539  else if(colldir1)
540  {
541  Blas::Dgemm('T','N',nmodes0,nquad1,nquad0,1.0,base0.get(),
542  nquad0,inarray.get(),nquad0,0.0,outarray.get(),nmodes0);
543  }
544  else
545  {
546  ASSERTL1(wsp.num_elements()>=nquad1*nmodes0,"Workspace size is not sufficient");
547 
548 #if 1
549  Blas::Dgemm('T','N',nmodes0,nquad1,nquad0,1.0,base0.get(),
550  nquad0,inarray.get(),nquad0,0.0,wsp.get(),nmodes0);
551 
552 #else
553  for(int i = 0; i < nmodes0; ++i)
554  {
555  for(int j = 0; j < nquad1; ++j)
556  {
557  wsp[j*nmodes0+i] = Blas::Ddot(nquad0,
558  base0.get()+i*nquad0,1,
559  inarray.get()+j*nquad0,1);
560  }
561  }
562 #endif
563  Blas::Dgemm('N', 'N',nmodes0,nmodes1, nquad1,1.0, wsp.get(),
564  nmodes0, base1.get(), nquad1, 0.0,outarray.get(),nmodes0);
565  }
566  }
567 
568  //////////////////////////
569  // Evaluation functions //
570  //////////////////////////
571 
572 
575  {
576  eta[0] = xi[0];
577  eta[1] = xi[1];
578  }
579 
580  /** \brief Fill outarray with mode \a mode of expansion
581  *
582  * Note for quadrilateral expansions _base[0] (i.e. p) modes run
583  * fastest
584  */
585 
586  void StdQuadExp::v_FillMode(const int mode,
587  Array<OneD, NekDouble> &outarray)
588  {
589  int i;
590  int nquad0 = m_base[0]->GetNumPoints();
591  int nquad1 = m_base[1]->GetNumPoints();
592  Array<OneD, const NekDouble> base0 = m_base[0]->GetBdata();
593  Array<OneD, const NekDouble> base1 = m_base[1]->GetBdata();
594  int btmp0 = m_base[0]->GetNumModes();
595  int mode0 = mode%btmp0;
596  int mode1 = mode/btmp0;
597 
598 
599  ASSERTL2(mode1 == (int)floor((1.0*mode)/btmp0),
600  "Integer Truncation not Equiv to Floor");
601 
602  ASSERTL2(m_ncoeffs <= mode,
603  "calling argument mode is larger than total expansion order");
604 
605  for(i = 0; i < nquad1; ++i)
606  {
607  Vmath::Vcopy(nquad0,(NekDouble *)(base0.get() + mode0*nquad0),
608  1, &outarray[0]+i*nquad0,1);
609  }
610 
611  for(i = 0; i < nquad0; ++i)
612  {
613  Vmath::Vmul(nquad1,(NekDouble *)(base1.get() + mode1*nquad1),1,
614  &outarray[0]+i,nquad0,&outarray[0]+i,nquad0);
615  }
616  }
617 
618  //////////////////////
619  // Helper functions //
620  //////////////////////
621 
623  {
624  return 4;
625  }
626 
628  {
629  return 4;
630  }
631 
632  int StdQuadExp::v_GetEdgeNcoeffs(const int i) const
633  {
634  ASSERTL2((i >= 0)&&(i <= 3),"edge id is out of range");
635 
636  if((i == 0)||(i == 2))
637  {
638  return GetBasisNumModes(0);
639  }
640  else
641  {
642  return GetBasisNumModes(1);
643  }
644  }
645 
646  int StdQuadExp::v_GetEdgeNumPoints(const int i) const
647  {
648  ASSERTL2((i >= 0)&&(i <= 3),"edge id is out of range");
649 
650  if((i == 0)||(i == 2))
651  {
652  return GetNumPoints(0);
653  }
654  else
655  {
656  return GetNumPoints(1);
657  }
658  }
659 
661  {
662  ASSERTL2((i >= 0)&&(i <= 3),"edge id is out of range");
663 
664  if((i == 0)||(i == 2))
665  {
666  return GetBasisType(0);
667  }
668  else
669  {
670  return GetBasisType(1);
671  }
672  }
673 
675  {
676  ASSERTL2((edge >= 0)&&(edge <= 3),"edge id is out of range");
677 
678  if((edge == 0)||(edge == 2))
679  {
680  return 0;
681  }
682  else
683  {
684  return 1;
685  }
686  }
687 
689  {
690  ASSERTL2((i >= 0)&&(i <= 3),"edge id is out of range");
691 
692  if((i == 0)||(i == 2))
693  {
694  return GetBasis(0)->GetBasisKey();
695  }
696  else
697  {
698  return GetBasis(1)->GetBasisKey();
699  }
700 
701  }
702 
704  {
706  }
707 
708 
710  {
714  "BasisType is not a boundary interior form");
718  "BasisType is not a boundary interior form");
719 
720  return 4 + 2*(GetBasisNumModes(0)-2) + 2*(GetBasisNumModes(1)-2);
721  }
722 
724  {
728  "BasisType is not a boundary interior form");
732  "BasisType is not a boundary interior form");
733 
734  return 2*GetBasisNumModes(0) + 2*GetBasisNumModes(1);
735  }
736 
738  const std::vector<unsigned int> &nummodes,
739  int &modes_offset)
740  {
741  int nmodes = nummodes[modes_offset]*nummodes[modes_offset+1];
742  modes_offset += 2;
743 
744  return nmodes;
745  }
746 
748  {
749  bool returnval = false;
750 
753  {
756  {
757  returnval = true;
758  }
759  }
760 
761  return returnval;
762  }
763 
765  Array<OneD, NekDouble> &coords_1,
766  Array<OneD, NekDouble> &coords_2)
767  {
768  boost::ignore_unused(coords_2);
769  Array<OneD, const NekDouble> z0 = m_base[0]->GetZ();
770  Array<OneD, const NekDouble> z1 = m_base[1]->GetZ();
771  int nq0 = GetNumPoints(0);
772  int nq1 = GetNumPoints(1);
773  int i;
774 
775  for(i = 0; i < nq1; ++i)
776  {
777  Blas::Dcopy(nq0,z0.get(), 1,&coords_0[0] + i*nq0,1);
778  Vmath::Fill(nq0,z1[i],&coords_1[0] + i*nq0,1);
779  }
780  }
781 
782  //////////////
783  // Mappings //
784  //////////////
785 
787  {
788  int i;
789  int cnt=0;
790  int nummodes0, nummodes1;
791  int value1 = 0, value2 = 0;
792  if(outarray.num_elements()!=NumBndryCoeffs())
793  {
795  }
796 
797  nummodes0 = m_base[0]->GetNumModes();
798  nummodes1 = m_base[1]->GetNumModes();
799 
800  const LibUtilities::BasisType Btype0 = GetBasisType(0);
801  const LibUtilities::BasisType Btype1 = GetBasisType(1);
802 
803  switch(Btype1)
804  {
807  value1 = nummodes0;
808  break;
810  value1 = 2*nummodes0;
811  break;
812  default:
813  ASSERTL0(0,"Mapping array is not defined for this expansion");
814  break;
815  }
816 
817  for(i = 0; i < value1; i++)
818  {
819  outarray[i]=i;
820  }
821  cnt=value1;
822 
823  switch(Btype0)
824  {
827  value2 = value1+nummodes0-1;
828  break;
830  value2 = value1+1;
831  break;
832  default:
833  ASSERTL0(0,"Mapping array is not defined for this expansion");
834  break;
835  }
836 
837  for(i = 0; i < nummodes1-2; i++)
838  {
839  outarray[cnt++]=value1+i*nummodes0;
840  outarray[cnt++]=value2+i*nummodes0;
841  }
842 
843 
845  {
846  for(i = nummodes0*(nummodes1-1);i < GetNcoeffs(); i++)
847  {
848  outarray[cnt++] = i;
849  }
850  }
851  }
852 
854  {
855  int i,j;
856  int cnt=0;
857  int nummodes0, nummodes1;
858  int startvalue = 0;
859  if(outarray.num_elements()!=GetNcoeffs()-NumBndryCoeffs())
860  {
862  }
863 
864  nummodes0 = m_base[0]->GetNumModes();
865  nummodes1 = m_base[1]->GetNumModes();
866 
867  const LibUtilities::BasisType Btype0 = GetBasisType(0);
868  const LibUtilities::BasisType Btype1 = GetBasisType(1);
869 
870  switch(Btype1)
871  {
873  startvalue = nummodes0;
874  break;
876  startvalue = 2*nummodes0;
877  break;
878  default:
879  ASSERTL0(0,"Mapping array is not defined for this expansion");
880  break;
881  }
882 
883  switch(Btype0)
884  {
886  startvalue++;
887  break;
889  startvalue+=2;
890  break;
891  default:
892  ASSERTL0(0,"Mapping array is not defined for this expansion");
893  break;
894  }
895 
896  for(i = 0; i < nummodes1-2; i++)
897  {
898  for(j = 0; j < nummodes0-2; j++)
899  {
900  outarray[cnt++]=startvalue+j;
901  }
902  startvalue+=nummodes0;
903  }
904  }
905 
906  int StdQuadExp::v_GetVertexMap(int localVertexId, bool useCoeffPacking)
907  {
908  int localDOF = 0;
909 
910  if(useCoeffPacking == true)
911  {
912  switch(localVertexId)
913  {
914  case 0:
915  {
916  localDOF = 0;
917  }
918  break;
919  case 1:
920  {
922  {
923  localDOF = m_base[0]->GetNumModes()-1;
924  }
925  else
926  {
927  localDOF = 1;
928  }
929  }
930  break;
931  case 2:
932  {
934  {
935  localDOF = m_base[0]->GetNumModes() * (m_base[1]->GetNumModes()-1);
936  }
937  else
938  {
939  localDOF = m_base[0]->GetNumModes();
940  }
941  }
942  break;
943  case 3:
944  {
946  {
947  localDOF = m_base[0]->GetNumModes()*m_base[1]->GetNumModes()-1;
948  }
949  else
950  {
951  localDOF = m_base[0]->GetNumModes()+1;
952  }
953  }
954  break;
955  default:
956  ASSERTL0(false,"eid must be between 0 and 3");
957  break;
958  }
959 
960  }
961  else
962  {
963  switch(localVertexId)
964  {
965  case 0:
966  {
967  localDOF = 0;
968  }
969  break;
970  case 1:
971  {
973  {
974  localDOF = m_base[0]->GetNumModes()-1;
975  }
976  else
977  {
978  localDOF = 1;
979  }
980  }
981  break;
982  case 2:
983  {
985  {
986  localDOF = m_base[0]->GetNumModes()*m_base[1]->GetNumModes()-1;
987  }
988  else
989  {
990  localDOF = m_base[0]->GetNumModes()+1;
991  }
992  }
993  break;
994  case 3:
995  {
997  {
998  localDOF = m_base[0]->GetNumModes() * (m_base[1]->GetNumModes()-1);
999  }
1000  else
1001  {
1002  localDOF = m_base[0]->GetNumModes();
1003  }
1004  }
1005  break;
1006  default:
1007  ASSERTL0(false,"eid must be between 0 and 3");
1008  break;
1009  }
1010  }
1011  return localDOF;
1012  }
1013 
1015  const Orientation edgeOrient,
1016  Array<OneD, unsigned int> &maparray,
1017  Array<OneD, int> &signarray)
1018  {
1019  int i;
1020  const int nummodes0 = m_base[0]->GetNumModes();
1021  const int nummodes1 = m_base[1]->GetNumModes();
1022  const int nEdgeIntCoeffs = GetEdgeNcoeffs(eid)-2;
1023  const LibUtilities::BasisType bType = GetEdgeBasisType(eid);
1024 
1025  if(maparray.num_elements() != nEdgeIntCoeffs)
1026  {
1027  maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
1028  }
1029 
1030  if(signarray.num_elements() != nEdgeIntCoeffs)
1031  {
1032  signarray = Array<OneD, int>(nEdgeIntCoeffs,1);
1033  }
1034  else
1035  {
1036  fill( signarray.get() , signarray.get()+nEdgeIntCoeffs, 1 );
1037  }
1038 
1039  if(bType == LibUtilities::eModified_A)
1040  {
1041  switch(eid)
1042  {
1043  case 0:
1044  {
1045  for(i = 0; i < nEdgeIntCoeffs; i++)
1046  {
1047  maparray[i] = i+2;
1048  }
1049 
1050  if(edgeOrient==eBackwards)
1051  {
1052  for(i = 1; i < nEdgeIntCoeffs; i+=2)
1053  {
1054  signarray[i] = -1;
1055  }
1056  }
1057  }
1058  break;
1059  case 1:
1060  {
1061  for(i = 0; i < nEdgeIntCoeffs; i++)
1062  {
1063  maparray[i] = (i+2)*nummodes0 + 1;
1064  }
1065 
1066  if(edgeOrient==eBackwards)
1067  {
1068  for(i = 1; i < nEdgeIntCoeffs; i+=2)
1069  {
1070  signarray[i] = -1;
1071  }
1072  }
1073  }
1074  break;
1075  case 2:
1076  {
1077  for(i = 0; i < nEdgeIntCoeffs; i++)
1078  {
1079  maparray[i] = nummodes0+i+2;
1080  }
1081 
1082  if(edgeOrient==eBackwards)
1083  {
1084  for(i = 1; i < nEdgeIntCoeffs; i+=2)
1085  {
1086  signarray[i] = -1;
1087  }
1088  }
1089  }
1090  break;
1091  case 3:
1092  {
1093  for(i = 0; i < nEdgeIntCoeffs; i++)
1094  {
1095  maparray[i] = (i+2)*nummodes0;
1096  }
1097 
1098  if(edgeOrient==eBackwards)
1099  {
1100  for(i = 1; i < nEdgeIntCoeffs; i+=2)
1101  {
1102  signarray[i] = -1;
1103  }
1104  }
1105  }
1106  break;
1107  default:
1108  ASSERTL0(false,"eid must be between 0 and 3");
1109  break;
1110  }
1111  }
1112  else if(bType == LibUtilities::eGLL_Lagrange)
1113  {
1114  switch(eid)
1115  {
1116  case 0:
1117  {
1118  for(i = 0; i < nEdgeIntCoeffs; i++)
1119  {
1120  maparray[i] = i+1;
1121  }
1122  }
1123  break;
1124  case 1:
1125  {
1126  for(i = 0; i < nEdgeIntCoeffs; i++)
1127  {
1128  maparray[i] = (i+2)*nummodes0 - 1;
1129  }
1130  }
1131  break;
1132  case 2:
1133  {
1134  for(i = 0; i < nEdgeIntCoeffs; i++)
1135  {
1136  maparray[i] = nummodes0*(nummodes1-1) + i + 1;
1137  }
1138  }
1139  break;
1140  case 3:
1141  {
1142  for(i = 0; i < nEdgeIntCoeffs; i++)
1143  {
1144  maparray[i] = nummodes0 * (i+1);
1145  }
1146  }
1147  break;
1148  default:
1149  ASSERTL0(false,"eid must be between 0 and 3");
1150  break;
1151  }
1152  if(edgeOrient == eBackwards)
1153  {
1154  reverse( maparray.get() , maparray.get()+nEdgeIntCoeffs );
1155  }
1156  }
1157  else
1158  {
1159  ASSERTL0(false,"Mapping not defined for this type of basis");
1160  }
1161 
1162  }
1163 
1165  const int eid,
1166  const Orientation edgeOrient,
1167  Array<OneD, unsigned int>& maparray,
1168  Array<OneD, int>& signarray,
1169  int P)
1170  {
1171  int i;
1172  int numModes=0;
1173  int order0 = m_base[0]->GetNumModes();
1174  int order1 = m_base[1]->GetNumModes();
1175 
1176  switch (eid)
1177  {
1178  case 0:
1179  case 2:
1180  numModes = order0;
1181  break;
1182  case 1:
1183  case 3:
1184  numModes = order1;
1185  break;
1186  default:
1187  ASSERTL0(false,"eid must be between 0 and 3");
1188  }
1189 
1190  bool checkForZeroedModes = false;
1191  if (P == -1)
1192  {
1193  P = numModes;
1194  }
1195  else if(P != numModes)
1196  {
1197  checkForZeroedModes = true;
1198  }
1199  const LibUtilities::BasisType bType = GetEdgeBasisType(eid);
1200 
1201 
1202  if (maparray.num_elements() != P)
1203  {
1204  maparray = Array<OneD, unsigned int>(P);
1205  }
1206 
1207  if(signarray.num_elements() != P)
1208  {
1209  signarray = Array<OneD, int>(P, 1);
1210  }
1211  else
1212  {
1213  fill(signarray.get(), signarray.get()+P, 1);
1214  }
1215 
1216  if (bType == LibUtilities::eModified_A)
1217  {
1218  switch (eid)
1219  {
1220  case 0:
1221  {
1222  for (i = 0; i < P; i++)
1223  {
1224  maparray[i] = i;
1225  }
1226 
1227  if (edgeOrient == eBackwards)
1228  {
1229  swap(maparray[0], maparray[1]);
1230 
1231  for(i = 3; i < P; i+=2)
1232  {
1233  signarray[i] = -1;
1234  }
1235  }
1236  }
1237  break;
1238  case 1:
1239  {
1240  for (i = 0; i < P; i++)
1241  {
1242  maparray[i] = i*order0 + 1;
1243  }
1244 
1245  if (edgeOrient == eBackwards)
1246  {
1247  swap(maparray[0], maparray[1]);
1248 
1249  for(i = 3; i < P; i+=2)
1250  {
1251  signarray[i] = -1;
1252  }
1253  }
1254  }
1255  break;
1256  case 2:
1257  {
1258  for (i = 0; i < P; i++)
1259  {
1260  maparray[i] = order0+i;
1261  }
1262 
1263  if (edgeOrient == eBackwards)
1264  {
1265  swap(maparray[0], maparray[1]);
1266 
1267  for (i = 3; i < P; i+=2)
1268  {
1269  signarray[i] = -1;
1270  }
1271  }
1272  }
1273  break;
1274  case 3:
1275  {
1276  for (i = 0; i < P; i++)
1277  {
1278  maparray[i] = i*order0;
1279  }
1280 
1281  if (edgeOrient == eBackwards)
1282  {
1283  swap(maparray[0], maparray[1]);
1284 
1285  for (i = 3; i < P; i+=2)
1286  {
1287  signarray[i] = -1;
1288  }
1289  }
1290  }
1291  break;
1292  default:
1293  ASSERTL0(false, "eid must be between 0 and 3");
1294  break;
1295  }
1296  }
1297  else if(bType == LibUtilities::eGLL_Lagrange ||
1299  {
1300  switch (eid)
1301  {
1302  case 0:
1303  {
1304  for (i = 0; i < P; i++)
1305  {
1306  maparray[i] = i;
1307  }
1308  }
1309  break;
1310  case 1:
1311  {
1312  for (i = 0; i < P; i++)
1313  {
1314  maparray[i] = (i+1)*order0 - 1;
1315  }
1316  }
1317  break;
1318  case 2:
1319  {
1320  for (i = 0; i < P; i++)
1321  {
1322  maparray[i] = order0*(order1-1) + i;
1323  }
1324  }
1325  break;
1326  case 3:
1327  {
1328  for (i = 0; i < P; i++)
1329  {
1330  maparray[i] = order0*i;
1331  }
1332  }
1333  break;
1334  default:
1335  ASSERTL0(false, "eid must be between 0 and 3");
1336  break;
1337  }
1338  if (edgeOrient == eBackwards)
1339  {
1340  reverse(maparray.get(), maparray.get()+P);
1341  }
1342  }
1343  else
1344  {
1345  ASSERTL0(false, "Mapping not defined for this type of basis");
1346  }
1347 
1348  if (checkForZeroedModes)
1349  {
1350  if (bType == LibUtilities::eModified_A)
1351  {
1352  // Zero signmap and set maparray to zero if
1353  // elemental modes are not as large as face modesl
1354  for (int j = numModes; j < P; j++)
1355  {
1356  signarray[j] = 0.0;
1357  maparray[j] = maparray[0];
1358  }
1359  }
1360  else
1361  {
1362  ASSERTL0(false, "Different trace space edge dimension "
1363  "and element edge dimension not possible "
1364  "for GLL-Lagrange bases");
1365  }
1366  }
1367  }
1368 
1369 
1370 
1371  ///////////////////////
1372  // Wrapper Functions //
1373  ///////////////////////
1374 
1376  {
1377  int i;
1378  int order0 = GetBasisNumModes(0);
1379  int order1 = GetBasisNumModes(1);
1380  MatrixType mtype = mkey.GetMatrixType();
1381 
1382  DNekMatSharedPtr Mat;
1383 
1384  switch(mtype)
1385  {
1387  {
1388  int nq0 = m_base[0]->GetNumPoints();
1389  int nq1 = m_base[1]->GetNumPoints();
1390  int nq;
1391 
1392  // take definition from key
1394  {
1395  nq = (int) mkey.GetConstFactor(eFactorConst);
1396  }
1397  else
1398  {
1399  nq = max(nq0,nq1);
1400  }
1401 
1402  int neq = LibUtilities::StdQuadData::
1403  getNumberOfCoefficients(nq, nq);
1404  Array<OneD, Array<OneD, NekDouble> > coords(neq);
1405  Array<OneD, NekDouble> coll (2);
1407  Array<OneD, NekDouble> tmp (nq0);
1408 
1410  AllocateSharedPtr(neq, nq0 * nq1);
1411  int cnt = 0;
1412 
1413  for(int i = 0; i < nq; ++i)
1414  {
1415  for(int j = 0; j < nq; ++j,++cnt)
1416  {
1417  coords[cnt] = Array<OneD, NekDouble>(2);
1418  coords[cnt][0] = -1.0 + 2*j/(NekDouble)(nq-1);
1419  coords[cnt][1] = -1.0 + 2*i/(NekDouble)(nq-1);
1420  }
1421  }
1422 
1423  for(int i = 0; i < neq; ++i)
1424  {
1425  LocCoordToLocCollapsed(coords[i],coll);
1426 
1427  I[0] = m_base[0]->GetI(coll);
1428  I[1] = m_base[1]->GetI(coll+1);
1429 
1430  // interpolate first coordinate direction
1431  for (int j = 0; j < nq1; ++j)
1432  {
1433  NekDouble fac = (I[1]->GetPtr())[j];
1434  Vmath::Smul(nq0,fac,I[0]->GetPtr(),1,tmp,1);
1435 
1436  Vmath::Vcopy(nq0, &tmp[0], 1,
1437  Mat->GetRawPtr()+j*nq0*neq+i,neq);
1438  }
1439 
1440  }
1441  break;
1442  }
1443  case eMass:
1444  {
1446  // For Fourier basis set the imaginary component of mean mode
1447  // to have a unit diagonal component in mass matrix
1449  {
1450  for(i = 0; i < order1; ++i)
1451  {
1452  (*Mat)(order0*i+1,i*order0+1) = 1.0;
1453  }
1454  }
1455 
1457  {
1458  for(i = 0; i < order0; ++i)
1459  {
1460  (*Mat)(order0+i ,order0+i) = 1.0;
1461  }
1462  }
1463  break;
1464  }
1465  case eFwdTrans:
1466  {
1468  StdMatrixKey iprodkey(eIProductWRTBase,DetShapeType(),*this);
1469  DNekMat &Iprod = *GetStdMatrix(iprodkey);
1470  StdMatrixKey imasskey(eInvMass,DetShapeType(),*this);
1471  DNekMat &Imass = *GetStdMatrix(imasskey);
1472 
1473  (*Mat) = Imass*Iprod;
1474  break;
1475  }
1476  case eGaussDG:
1477  {
1478  ConstFactorMap factors = mkey.GetConstFactors();
1479 
1480  int edge = (int)factors[StdRegions::eFactorGaussEdge];
1481  int dir = (edge + 1) % 2;
1482  int nCoeffs = m_base[dir]->GetNumModes();
1483 
1484  const LibUtilities::PointsKey BS_p(
1486  const LibUtilities::BasisKey BS_k(
1487  LibUtilities::eGauss_Lagrange, nCoeffs, BS_p);
1488 
1489  Array<OneD, NekDouble> coords(1, 0.0);
1490  coords[0] = (edge == 0 || edge == 3) ? -1.0 : 1.0;
1491 
1494  DNekMatSharedPtr m_Ix = basis->GetI(coords);
1495 
1497  1.0, nCoeffs);
1498  Vmath::Vcopy(nCoeffs, m_Ix->GetPtr(), 1, Mat->GetPtr(), 1);
1499  break;
1500  }
1501  default:
1502  {
1504  break;
1505  }
1506  }
1507 
1508  return Mat;
1509  }
1510 
1512  {
1513  return GenMatrix(mkey);
1514  }
1515 
1516  ///////////////////////////////////
1517  // Operator evaluation functions //
1518  ///////////////////////////////////
1519 
1521  const Array<OneD, const NekDouble> &inarray,
1522  Array<OneD,NekDouble> &outarray,
1523  const StdMatrixKey &mkey)
1524  {
1525 
1527 
1528  if(inarray.get() == outarray.get())
1529  {
1531  Vmath::Vcopy(m_ncoeffs,inarray.get(),1,tmp.get(),1);
1532 
1533  Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, 1.0, mat->GetPtr().get(),
1534  m_ncoeffs, tmp.get(), 1, 0.0, outarray.get(), 1);
1535  }
1536  else
1537  {
1538  Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, 1.0, mat->GetPtr().get(),
1539  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
1540  }
1541  }
1542 
1543 
1545  const StdMatrixKey &mkey)
1546  {
1547  // Generate an orthonogal expansion
1548  int qa = m_base[0]->GetNumPoints();
1549  int qb = m_base[1]->GetNumPoints();
1550  int nmodes_a = m_base[0]->GetNumModes();
1551  int nmodes_b = m_base[1]->GetNumModes();
1552  int nmodes = min(nmodes_a,nmodes_b);
1553  // Declare orthogonal basis.
1556 
1559  StdQuadExp OrthoExp(Ba,Bb);
1560 
1561  //SVV parameters loaded from the .xml case file
1562  Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
1563 
1564  // project onto modal space.
1565  OrthoExp.FwdTrans(array,orthocoeffs);
1566 
1567  if(mkey.ConstFactorExists(eFactorSVVPowerKerDiffCoeff)) // Rodrigo's power kernel
1568  {
1570  NekDouble SvvDiffCoeff =
1573 
1574  for(int j = 0; j < nmodes_a; ++j)
1575  {
1576  for(int k = 0; k < nmodes_b; ++k)
1577  {
1578  // linear space but makes high modes very negative
1579  NekDouble fac = std::max(
1580  pow((1.0*j)/(nmodes_a-1),cutoff*nmodes_a),
1581  pow((1.0*k)/(nmodes_b-1),cutoff*nmodes_b));
1582 
1583  orthocoeffs[j*nmodes_b+k] *= SvvDiffCoeff * fac;
1584  }
1585  }
1586  }
1587  else if(mkey.ConstFactorExists(eFactorSVVDGKerDiffCoeff)) // Rodrigo/mansoor's DG kernel
1588  {
1589  NekDouble SvvDiffCoeff =
1592  int max_ab = max(nmodes_a-kSVVDGFiltermodesmin,
1593  nmodes_b-kSVVDGFiltermodesmin);
1594  max_ab = max(max_ab,0);
1595  max_ab = min(max_ab,kSVVDGFiltermodesmax-kSVVDGFiltermodesmin);
1596 
1597  for(int j = 0; j < nmodes_a; ++j)
1598  {
1599  for(int k = 0; k < nmodes_b; ++k)
1600  {
1601  int maxjk = max(j,k);
1602  maxjk = min(maxjk,kSVVDGFiltermodesmax-1);
1603 
1604  orthocoeffs[j*nmodes_b+k] *= SvvDiffCoeff *
1605  kSVVDGFilter[max_ab][maxjk];
1606  }
1607  }
1608  }
1609  else
1610  {
1611  NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDiffCoeff);
1612  //Exponential Kernel implementation
1613  int cutoff = (int) (mkey.GetConstFactor(eFactorSVVCutoffRatio)*
1614  min(nmodes_a,nmodes_b));
1615 
1616  //counters for scanning through orthocoeffs array
1617  int cnt = 0;
1618 
1619  //------"New" Version August 22nd '13--------------------
1620  for(int j = 0; j < nmodes_a; ++j)
1621  {
1622  for(int k = 0; k < nmodes_b; ++k)
1623  {
1624  if(j + k >= cutoff)//to filter out only the "high-modes"
1625  {
1626  orthocoeffs[j*nmodes_b+k] *=
1627  (SvvDiffCoeff*exp(-(j+k-nmodes)*(j+k-nmodes)/
1628  ((NekDouble)((j+k-cutoff+1)*
1629  (j+k-cutoff+1)))));
1630  }
1631  else
1632  {
1633  orthocoeffs[j*nmodes_b+k] *= 0.0;
1634  }
1635  cnt++;
1636  }
1637  }
1638  }
1639 
1640  // backward transform to physical space
1641  OrthoExp.BwdTrans(orthocoeffs,array);
1642  }
1643 
1645  Array<OneD, NekDouble> &array,
1646  const NekDouble alpha,
1647  const NekDouble exponent,
1648  const NekDouble cutoff)
1649  {
1650  // Generate an orthogonal expansion
1651  int qa = m_base[0]->GetNumPoints();
1652  int qb = m_base[1]->GetNumPoints();
1653  int nmodesA = m_base[0]->GetNumModes();
1654  int nmodesB = m_base[1]->GetNumModes();
1655  int P = nmodesA - 1;
1656  int Q = nmodesB - 1;
1657 
1658  // Declare orthogonal basis.
1661 
1664  StdQuadExp OrthoExp(Ba,Bb);
1665 
1666  // Cutoff
1667  int Pcut = cutoff*P;
1668  int Qcut = cutoff*Q;
1669 
1670  // Project onto orthogonal space.
1671  Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
1672  OrthoExp.FwdTrans(array,orthocoeffs);
1673 
1674  //
1675  NekDouble fac, fac1, fac2;
1676  for(int i = 0; i < nmodesA; ++i)
1677  {
1678  for(int j = 0; j < nmodesB; ++j)
1679  {
1680  //to filter out only the "high-modes"
1681  if(i > Pcut || j > Qcut)
1682  {
1683  fac1 = (NekDouble) (i - Pcut)/( (NekDouble)(P - Pcut) );
1684  fac2 = (NekDouble) (j - Qcut)/( (NekDouble)(Q - Qcut) );
1685  fac = max(fac1, fac2);
1686  fac = pow(fac, exponent);
1687  orthocoeffs[i*nmodesB+j] *= exp(-alpha*fac);
1688  }
1689  }
1690  }
1691 
1692  // backward transform to physical space
1693  OrthoExp.BwdTrans(orthocoeffs,array);
1694  }
1695 
1697  int numMin,
1698  const Array<OneD, const NekDouble> &inarray,
1699  Array<OneD, NekDouble> &outarray)
1700  {
1701  int n_coeffs = inarray.num_elements();
1702 
1703 
1704  Array<OneD, NekDouble> coeff(n_coeffs);
1705  Array<OneD, NekDouble> coeff_tmp(n_coeffs,0.0);
1708 
1709  int nmodes0 = m_base[0]->GetNumModes();
1710  int nmodes1 = m_base[1]->GetNumModes();
1711  int numMax = nmodes0;
1712 
1713  Vmath::Vcopy(n_coeffs,inarray,1,coeff_tmp,1);
1714 
1715  const LibUtilities::PointsKey Pkey0(
1717  const LibUtilities::PointsKey Pkey1(
1719 
1720  LibUtilities::BasisKey b0(m_base[0]->GetBasisType(),nmodes0,Pkey0);
1721  LibUtilities::BasisKey b1(m_base[1]->GetBasisType(),nmodes1,Pkey1);
1722 
1723  LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A,nmodes0,Pkey0);
1724  LibUtilities::BasisKey bortho1(LibUtilities::eOrtho_A,nmodes1,Pkey1);
1725 
1727  b0, b1, coeff_tmp, bortho0, bortho1, coeff);
1728 
1729  Vmath::Zero(n_coeffs,coeff_tmp,1);
1730 
1731  int cnt = 0;
1732  for (int i = 0; i < numMin+1; ++i)
1733  {
1734  Vmath::Vcopy(numMin,
1735  tmp = coeff+cnt,1,
1736  tmp2 = coeff_tmp+cnt,1);
1737 
1738  cnt = i*numMax;
1739  }
1740 
1742  bortho0, bortho1, coeff_tmp, b0, b1, outarray);
1743  }
1744 
1745 
1747  const Array<OneD, const NekDouble> &inarray,
1748  Array<OneD,NekDouble> &outarray,
1749  const StdMatrixKey &mkey)
1750  {
1751  StdExpansion::MassMatrixOp_MatFree(inarray,outarray,mkey);
1752  }
1753 
1755  const Array<OneD, const NekDouble> &inarray,
1756  Array<OneD,NekDouble> &outarray,
1757  const StdMatrixKey &mkey)
1758  {
1759  StdQuadExp::v_LaplacianMatrixOp_MatFree(inarray,outarray,mkey);
1760  }
1761 
1762  void StdQuadExp::v_LaplacianMatrixOp(const int k1, const int k2,
1763  const Array<OneD, const NekDouble> &inarray,
1764  Array<OneD,NekDouble> &outarray,
1765  const StdMatrixKey &mkey)
1766  {
1767  StdExpansion::LaplacianMatrixOp_MatFree(k1,k2,inarray,outarray,mkey);
1768  }
1769 
1771  const Array<OneD, const NekDouble> &inarray,
1772  Array<OneD,NekDouble> &outarray,
1773  const StdMatrixKey &mkey)
1774  {
1775  StdExpansion::WeakDerivMatrixOp_MatFree(i,inarray,outarray,mkey);
1776  }
1777 
1779  Array<OneD,NekDouble> &outarray,
1780  const StdMatrixKey &mkey)
1781  {
1782  StdQuadExp::v_HelmholtzMatrixOp_MatFree(inarray,outarray,mkey);
1783  }
1784 
1785  //up to here
1787  const Array<OneD, const NekDouble> &inarray,
1788  Array<OneD, NekDouble> &outarray)
1789  {
1790  int i;
1791  int nquad0 = m_base[0]->GetNumPoints();
1792  int nquad1 = m_base[1]->GetNumPoints();
1793 
1794  const Array<OneD, const NekDouble>& w0 = m_base[0]->GetW();
1795  const Array<OneD, const NekDouble>& w1 = m_base[1]->GetW();
1796 
1797  // multiply by integration constants
1798  for(i = 0; i < nquad1; ++i)
1799  {
1800  Vmath::Vmul(nquad0, inarray.get()+i*nquad0,1,
1801  w0.get(),1,outarray.get()+i*nquad0,1);
1802  }
1803 
1804  for(i = 0; i < nquad0; ++i)
1805  {
1806  Vmath::Vmul(nquad1,outarray.get()+i,nquad0,w1.get(),1,
1807  outarray.get()+i,nquad0);
1808  }
1809  }
1810 
1812  Array<OneD, int> &conn,
1813  bool standard)
1814  {
1815  boost::ignore_unused(standard);
1816 
1817  int np1 = m_base[0]->GetNumPoints();
1818  int np2 = m_base[1]->GetNumPoints();
1819  int np = max(np1,np2);
1820 
1821  conn = Array<OneD, int>(6*(np-1)*(np-1));
1822 
1823  int row = 0;
1824  int rowp1 = 0;
1825  int cnt = 0;
1826  for(int i = 0; i < np-1; ++i)
1827  {
1828  rowp1 += np;
1829  for(int j = 0; j < np-1; ++j)
1830  {
1831  conn[cnt++] = row +j;
1832  conn[cnt++] = row +j+1;
1833  conn[cnt++] = rowp1 +j;
1834 
1835  conn[cnt++] = rowp1 +j+1;
1836  conn[cnt++] = rowp1 +j;
1837  conn[cnt++] = row +j+1;
1838  }
1839  row += np;
1840  }
1841  }
1842 
1843  } //end namespace
1844 }//end namespace
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey)
virtual int v_GetNverts() const
Definition: StdQuadExp.cpp:622
DNekMatSharedPtr GenMatrix(const StdMatrixKey &mkey)
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
virtual LibUtilities::BasisType v_GetEdgeBasisType(const int i) const
Definition: StdQuadExp.cpp:660
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdQuadExp.cpp:447
const int kSVVDGFiltermodesmax
Definition: StdRegions.hpp:385
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:974
static Array< OneD, NekDouble > NullNekDouble1DArray
void BwdTrans_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:469
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:16
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Calculate the derivative of the physical points.
Definition: StdQuadExp.cpp:99
void IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:228
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:945
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey)
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:177
void IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
this function calculates the inner product of a given function f with the different modes of the expa...
Definition: StdExpansion.h:634
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45
virtual void v_ReduceOrderCoeffs(int numMin, 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)
Principle Modified Functions .
Definition: BasisType.h:48
virtual bool v_IsBoundaryInteriorExpansion()
Definition: StdQuadExp.cpp:747
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:55
STL namespace.
virtual int v_GetEdgeNumPoints(const int i) const
Definition: StdQuadExp.cpp:646
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:294
NekDouble Integral(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &w0, const Array< OneD, const NekDouble > &w1)
virtual void v_GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdQuadExp.cpp:853
virtual void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
Fourier Expansion .
Definition: BasisType.h:53
const int kSVVDGFiltermodesmin
Definition: StdRegions.hpp:384
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
LibUtilities::NekManager< StdMatrixKey, DNekBlkMat, StdMatrixKey::opLess > m_stdStaticCondMatrixManager
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:81
std::shared_ptr< Basis > BasisSharedPtr
static void Dgemm(const char &transa, const char &transb, const int &m, const int &n, const int &k, const double &alpha, const double *a, const int &lda, const double *b, const int &ldb, const double &beta, double *c, const int &ldc)
BLAS level 3: Matrix-matrix multiply C = A x B where A[m x n], B[n x k], C[m x k].
Definition: Blas.hpp:213
virtual void v_StdPhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Definition: StdQuadExp.cpp:132
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:714
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d0, Array< OneD, NekDouble > &outarray_d1)
Calculate the 2D derivative in the local tensor/collapsed coordinate at the physical points...
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:48
virtual void v_FwdTrans_BndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdQuadExp.cpp:260
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Transform a given function from physical quadrature space to coefficient space.
Definition: StdQuadExp.cpp:237
void v_GeneralMatrixOp_MatOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
BasisManagerT & BasisManager(void)
virtual void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray)
Definition: StdQuadExp.cpp:786
virtual LibUtilities::ShapeType v_DetShapeType() const
Definition: StdQuadExp.cpp:703
void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1)
Definition: StdQuadExp.cpp:513
The base class for all shapes.
Definition: StdExpansion.h:68
int GetEdgeNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th edge.
Definition: StdExpansion.h:286
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:216
virtual void v_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
const LibUtilities::BasisSharedPtr & GetBasis(int dir) const
This function gets the shared point to basis in the dir direction.
Definition: StdExpansion.h:117
virtual void v_IProductWRTBase_MatOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdQuadExp.cpp:428
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdQuadExp.cpp:170
virtual int v_DetCartesianDirOfEdge(const int edge)
Definition: StdQuadExp.cpp:674
virtual void v_ExponentialFilter(Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
void WeakDerivMatrixOp_MatFree(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:817
virtual void v_GetEdgeInteriorMap(const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
Principle Orthogonal Functions .
Definition: BasisType.h:45
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdQuadExp.cpp:440
void v_GetEdgeToElementMap(const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int P=-1)
virtual void v_LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Definition: StdQuadExp.cpp:573
Defines a specification for a set of points.
Definition: Points.h:59
virtual void v_WeakDerivMatrixOp(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
double NekDouble
static void Dgemv(const char &trans, const int &m, const int &n, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = A x where A[m x n].
Definition: Blas.hpp:168
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_0, Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2)
Definition: StdQuadExp.cpp:764
NekDouble GetConstFactor(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:121
virtual void v_BwdTrans_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1)
Definition: StdQuadExp.cpp:192
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:130
virtual int v_NumBndryCoeffs() const
Definition: StdQuadExp.cpp:709
const NekDouble kSVVDGFilter[9][11]
Definition: StdRegions.hpp:387
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:346
virtual void v_IProductWRTDerivBase_MatOp(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdQuadExp.cpp:478
const ConstFactorMap & GetConstFactors() const
Definition: StdMatrixKey.h:135
static double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = .
Definition: Blas.hpp:140
virtual int v_CalcNumberOfCoefficients(const std::vector< unsigned int > &nummodes, int &modes_offset)
Definition: StdQuadExp.cpp:737
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:215
void GetEdgeToElementMap(const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int P=-1)
Definition: StdExpansion.h:849
bool ConstFactorExists(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:130
int getNumberOfCoefficients(int Na, int Nb)
Definition: ShapeType.hpp:137
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:274
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:164
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Calculate the inner product of inarray with respect to the basis B=base0*base1 and put into outarray...
Definition: StdQuadExp.cpp:385
void MassMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual int v_GetVertexMap(int localVertexId, bool useCoeffPacking=false)
Definition: StdQuadExp.cpp:906
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Backward transformation from coefficient space to physical space...
Definition: StdExpansion.h:530
virtual int v_GetNedges() const
Definition: StdQuadExp.cpp:627
virtual const LibUtilities::BasisKey v_DetEdgeBasisKey(const int i) const
Definition: StdQuadExp.cpp:688
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
Definition: StdQuadExp.cpp:399
LibUtilities::NekManager< StdMatrixKey, DNekMat, StdMatrixKey::opLess > m_stdMatrixManager
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs
NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray)
Integrates the specified function over the domain.
Definition: StdQuadExp.cpp:80
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdQuadExp.cpp:156
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140
virtual int v_GetEdgeNcoeffs(const int i) const
Definition: StdQuadExp.cpp:632
virtual void v_FillMode(const int mode, Array< OneD, NekDouble > &array)
Fill outarray with mode mode of expansion.
Definition: StdQuadExp.cpp:586
void IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
virtual void v_MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Lagrange for SEM basis .
Definition: BasisType.h:54
std::shared_ptr< StdSegExp > StdSegExpSharedPtr
Definition: StdSegExp.h:46
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:376
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void InterpCoeff2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
Definition: InterpCoeff.cpp:72
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
static void Dcopy(const int &n, const double *x, const int &incx, double *y, const int &incy)
BLAS level 1: Copy x to y.
Definition: Blas.hpp:103
LibUtilities::BasisType GetEdgeBasisType(const int i) const
This function returns the type of expansion basis on the i-th edge.
Definition: StdExpansion.h:412
Describes the specification for a Basis.
Definition: Basis.h:49
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:186
virtual void v_GetSimplexEquiSpacedConnectivity(Array< OneD, int > &conn, bool standard=true)
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
virtual int v_NumDGBndryCoeffs() const
Definition: StdQuadExp.cpp:723
virtual void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Forward transformation from physical space to coefficient space...