Nektar++
StdTriExp.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File StdTriExp.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: Triangle routines built upon StdExpansion2D
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <boost/core/ignore_unused.hpp>
36 
38 #include <StdRegions/StdTriExp.h>
40 #include <StdRegions/StdSegExp.h> // for StdSegExp, etc
41 
42 using namespace std;
43 
44 namespace Nektar
45 {
46  namespace StdRegions
47  {
48 
49  StdTriExp::StdTriExp()
50  {
51  }
52 
53 
54  StdTriExp::StdTriExp(
55  const LibUtilities::BasisKey &Ba,
56  const LibUtilities::BasisKey &Bb) :
57  StdExpansion (LibUtilities::StdTriData::getNumberOfCoefficients(
58  Ba.GetNumModes(),
59  Bb.GetNumModes()),
60  2,Ba,Bb),
61  StdExpansion2D(LibUtilities::StdTriData::getNumberOfCoefficients(
62  Ba.GetNumModes(),
63  Bb.GetNumModes()),
64  Ba,Bb)
65  {
66  ASSERTL0(Ba.GetNumModes() <= Bb.GetNumModes(),
67  "order in 'a' direction is higher than order "
68  "in 'b' direction");
69  }
70 
72  StdExpansion(T),
74  {
75  }
76 
78  {
79  }
80 
81  //-------------------------------
82  // Integration Methods
83  //-------------------------------
85  const Array<OneD, const NekDouble>& inarray)
86  {
87  int i;
88  int nquad1 = m_base[1]->GetNumPoints();
89  Array<OneD, NekDouble> w1_tmp(nquad1);
90 
91  Array<OneD, const NekDouble> w0 = m_base[0]->GetW();
92  Array<OneD, const NekDouble> w1 = m_base[1]->GetW();
93  Array<OneD, const NekDouble> z1 = m_base[1]->GetZ();
94 
95  switch(m_base[1]->GetPointsType())
96  {
97  case LibUtilities::eGaussRadauMAlpha1Beta0: // (0,1) Jacobi Inner product
98  {
99  Vmath::Smul(nquad1, 0.5, w1, 1, w1_tmp,1);
100  break;
101  }
102  default:
103  {
104  // include jacobian factor on whatever coordinates are defined.
105  for(i = 0; i < nquad1; ++i)
106  {
107  w1_tmp[i] = 0.5*(1-z1[i])*w1[i];
108  }
109  break;
110  }
111  }
112 
113  return StdExpansion2D::Integral(inarray,w0,w1_tmp);
114  }
115 
116  //-----------------------------
117  // Differentiation Methods
118  //-----------------------------
119 
120  /**
121  * \brief Calculate the derivative of the physical points.
122  *
123  * \f$ \frac{\partial u}{\partial x_1} = \left .
124  * \frac{2.0}{1-\eta_2} \frac{\partial u}{\partial d\eta_1}
125  * \right |_{\eta_2}\f$
126  *
127  * \f$ \frac{\partial u}{\partial x_2} = \left .
128  * \frac{1+\eta_1}{1-\eta_2} \frac{\partial u}{\partial d\eta_1}
129  * \right |_{\eta_2} + \left . \frac{\partial u}{\partial d\eta_2}
130  * \right |_{\eta_1} \f$
131  */
133  const Array<OneD, const NekDouble>& inarray,
134  Array<OneD, NekDouble>& out_d0,
135  Array<OneD, NekDouble>& out_d1,
136  Array<OneD, NekDouble>& out_d2)
137  {
138  boost::ignore_unused(out_d2);
139 
140  int i;
141  int nquad0 = m_base[0]->GetNumPoints();
142  int nquad1 = m_base[1]->GetNumPoints();
143  Array<OneD, NekDouble> wsp(std::max(nquad0, nquad1));
144 
145  const Array<OneD, const NekDouble>& z0 = m_base[0]->GetZ();
146  const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
147 
148  // set up geometric factor: 2/(1-z1)
149  Vmath::Sadd(nquad1, -1.0, z1, 1, wsp, 1);
150  Vmath::Sdiv(nquad1, -2.0, wsp, 1, wsp, 1);
151 
152  if (out_d0.num_elements() > 0)
153  {
154  PhysTensorDeriv(inarray, out_d0, out_d1);
155 
156  for (i = 0; i < nquad1; ++i)
157  {
158  Blas::Dscal(nquad0,wsp[i],&out_d0[0]+i*nquad0,1);
159  }
160 
161  // if no d1 required do not need to calculate both deriv
162  if (out_d1.num_elements() > 0)
163  {
164  // set up geometric factor: (1_z0)/(1-z1)
165  Vmath::Sadd(nquad0, 1.0, z0, 1, wsp, 1);
166  Vmath::Smul(nquad0, 0.5, wsp, 1, wsp, 1);
167 
168  for (i = 0; i < nquad1; ++i)
169  {
170  Vmath::Vvtvp(nquad0,&wsp[0],1,&out_d0[0]+i*nquad0,
171  1,&out_d1[0]+i*nquad0,
172  1,&out_d1[0]+i*nquad0,1);
173  }
174  }
175  }
176  else if (out_d1.num_elements() > 0)
177  {
178  Array<OneD, NekDouble> diff0(nquad0*nquad1);
179  PhysTensorDeriv(inarray, diff0, out_d1);
180 
181  for (i = 0; i < nquad1; ++i)
182  {
183  Blas::Dscal(nquad0,wsp[i],&diff0[0]+i*nquad0,1);
184  }
185 
186  Vmath::Sadd(nquad0, 1.0, z0, 1, wsp, 1);
187  Vmath::Smul(nquad0, 0.5, wsp, 1, wsp, 1);
188 
189  for (i = 0; i < nquad1; ++i)
190  {
191  Vmath::Vvtvp(nquad0,&wsp[0],1,&diff0[0]+i*nquad0,
192  1,&out_d1[0]+i*nquad0,
193  1,&out_d1[0]+i*nquad0,1);
194  }
195  }
196  }
197 
199  const int dir,
200  const Array<OneD, const NekDouble>& inarray,
201  Array<OneD, NekDouble>& outarray)
202  {
203  switch(dir)
204  {
205  case 0:
206  {
207  v_PhysDeriv(inarray, outarray, NullNekDouble1DArray);
208  break;
209  }
210  case 1:
211  {
212  v_PhysDeriv(inarray, NullNekDouble1DArray, outarray);
213  break;
214  }
215  default:
216  {
217  ASSERTL1(false,"input dir is out of range");
218  break;
219  }
220  }
221  }
222 
224  const Array<OneD, const NekDouble>& inarray,
225  Array<OneD, NekDouble>& out_d0,
226  Array<OneD, NekDouble>& out_d1,
227  Array<OneD, NekDouble>& out_d2)
228  {
229  boost::ignore_unused(out_d2);
230  StdTriExp::v_PhysDeriv(inarray, out_d0, out_d1);
231  }
232 
234  const int dir,
235  const Array<OneD, const NekDouble>& inarray,
236  Array<OneD, NekDouble>& outarray)
237  {
238  StdTriExp::v_PhysDeriv(dir,inarray,outarray);
239  }
240 
241 
242  //---------------------------------------
243  // Transforms
244  //---------------------------------------
245 
246  /**
247  * \brief Backward tranform for triangular elements
248  *
249  * @note 'q' (base[1]) runs fastest in this element.
250  */
252  const Array<OneD, const NekDouble>& inarray,
253  Array<OneD, NekDouble>& outarray)
254  {
255  v_BwdTrans_SumFac(inarray,outarray);
256  }
257 
258 
260  const Array<OneD, const NekDouble>& inarray,
261  Array<OneD, NekDouble>& outarray)
262  {
264  m_base[1]->GetNumModes());
265 
266  BwdTrans_SumFacKernel(m_base[0]->GetBdata(),
267  m_base[1]->GetBdata(),
268  inarray,outarray,wsp);
269  }
270 
272  const Array<OneD, const NekDouble>& base0,
273  const Array<OneD, const NekDouble>& base1,
274  const Array<OneD, const NekDouble>& inarray,
275  Array<OneD, NekDouble>& outarray,
277  bool doCheckCollDir0,
278  bool doCheckCollDir1)
279  {
280  boost::ignore_unused(doCheckCollDir0, doCheckCollDir1);
281 
282  int i;
283  int mode;
284  int nquad0 = m_base[0]->GetNumPoints();
285  int nquad1 = m_base[1]->GetNumPoints();
286  int nmodes0 = m_base[0]->GetNumModes();
287  int nmodes1 = m_base[1]->GetNumModes();
288 
289  ASSERTL1(wsp.num_elements() >= nquad0*nmodes1,
290  "Workspace size is not sufficient");
293  "Basis[1] is not of general tensor type");
294 
295  for (i = mode = 0; i < nmodes0; ++i)
296  {
297  Blas::Dgemv('N', nquad1,nmodes1-i,1.0,base1.get()+mode*nquad1,
298  nquad1,&inarray[0]+mode,1,0.0,&wsp[0]+i*nquad1,1);
299  mode += nmodes1-i;
300  }
301 
302  // fix for modified basis by splitting top vertex mode
304  {
305  Blas::Daxpy(nquad1,inarray[1],base1.get()+nquad1,1,
306  &wsp[0]+nquad1,1);
307  }
308 
309  Blas::Dgemm('N','T', nquad0,nquad1,nmodes0,1.0, base0.get(),nquad0,
310  &wsp[0], nquad1,0.0, &outarray[0], nquad0);
311  }
312 
314  const Array<OneD, const NekDouble>& inarray,
315  Array<OneD, NekDouble>& outarray)
316  {
317  v_IProductWRTBase(inarray,outarray);
318 
319  // get Mass matrix inverse
320  StdMatrixKey masskey(eInvMass,DetShapeType(),*this);
321  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
322 
323  // copy inarray in case inarray == outarray
324  NekVector<NekDouble> in(m_ncoeffs,outarray,eCopy);
326 
327  out = (*matsys)*in;
328  }
329 
330 
332  const Array<OneD, const NekDouble>& inarray,
333  Array<OneD, NekDouble>& outarray)
334  {
335  int i,j;
336  int npoints[2] = {m_base[0]->GetNumPoints(),
337  m_base[1]->GetNumPoints()};
338  int nmodes[2] = {m_base[0]->GetNumModes(),
339  m_base[1]->GetNumModes()};
340 
341  fill(outarray.get(), outarray.get()+m_ncoeffs, 0.0 );
342 
343  Array<OneD, NekDouble> physEdge[3];
344  Array<OneD, NekDouble> coeffEdge[3];
345  for(i = 0; i < 3; i++)
346  {
347  physEdge[i] = Array<OneD, NekDouble>(npoints[i!=0]);
348  coeffEdge[i] = Array<OneD, NekDouble>(nmodes[i!=0]);
349  }
350 
351  for(i = 0; i < npoints[0]; i++)
352  {
353  physEdge[0][i] = inarray[i];
354  }
355 
356  for(i = 0; i < npoints[1]; i++)
357  {
358  physEdge[1][i] = inarray[npoints[0]-1+i*npoints[0]];
359  physEdge[2][i] = inarray[(npoints[1]-1)*npoints[0]-i*npoints[0]];
360  }
361 
362  StdSegExpSharedPtr segexp[2] = {
364  m_base[0]->GetBasisKey()),
366  m_base[1]->GetBasisKey())
367  };
368 
369  Array<OneD, unsigned int> mapArray;
370  Array<OneD, int> signArray;
371  NekDouble sign;
372 
373  for (i = 0; i < 3; i++)
374  {
375  //segexp[i!=0]->v_FwdTrans_BndConstrained(physEdge[i],coeffEdge[i]);
376  segexp[i!=0]->FwdTrans_BndConstrained(physEdge[i],coeffEdge[i]);
377 
378  v_GetEdgeToElementMap(i,eForwards,mapArray,signArray);
379  for (j = 0; j < nmodes[i != 0]; j++)
380  {
381  sign = (NekDouble) signArray[j];
382  outarray[ mapArray[j] ] = sign * coeffEdge[i][j];
383  }
384  }
385 
388 
389  StdMatrixKey masskey(eMass,DetShapeType(),*this);
390  MassMatrixOp(outarray,tmp0,masskey);
391  v_IProductWRTBase(inarray,tmp1);
392 
393  Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
394 
395  // get Mass matrix inverse (only of interior DOF)
396  // use block (1,1) of the static condensed system
397  // note: this block alreay contains the inverse matrix
398  DNekMatSharedPtr matsys =
399  (m_stdStaticCondMatrixManager[masskey])->GetBlock(1,1);
400 
401  int nBoundaryDofs = v_NumBndryCoeffs();
402  int nInteriorDofs = m_ncoeffs - nBoundaryDofs;
403 
404  Array<OneD, NekDouble> rhs (nInteriorDofs);
405  Array<OneD, NekDouble> result(nInteriorDofs);
406 
407  v_GetInteriorMap(mapArray);
408 
409  for (i = 0; i < nInteriorDofs; i++)
410  {
411  rhs[i] = tmp1[ mapArray[i] ];
412  }
413 
414  Blas::Dgemv('N',nInteriorDofs,nInteriorDofs,
415  1.0,&(matsys->GetPtr())[0],nInteriorDofs,
416  rhs.get(),1,
417  0.0,result.get(),1);
418 
419  for (i = 0; i < nInteriorDofs; i++)
420  {
421  outarray[ mapArray[i] ] = result[i];
422  }
423  }
424 
425  //---------------------------------------
426  // Inner product functions
427  //---------------------------------------
428 
429  /**
430  * \brief Calculate the inner product of inarray with respect to the
431  * basis B=base0[p]*base1[pq] and put into outarray.
432  *
433  * \f$
434  * \begin{array}{rcl}
435  * I_{pq} = (\phi^A_q \phi^B_{pq}, u) &=&
436  * \sum_{i=0}^{nq_0}\sum_{j=0}^{nq_1}
437  * \phi^A_p(\eta_{0,i})\phi^B_{pq}(\eta_{1,j}) w^0_i w^1_j
438  * u(\xi_{0,i} \xi_{1,j}) \\
439  * & = & \sum_{i=0}^{nq_0} \phi^A_p(\eta_{0,i})
440  * \sum_{j=0}^{nq_1} \phi^B_{pq}(\eta_{1,j}) \tilde{u}_{i,j}
441  * \end{array}
442  * \f$
443  *
444  * where
445  *
446  * \f$ \tilde{u}_{i,j} = w^0_i w^1_j u(\xi_{0,i},\xi_{1,j}) \f$
447  *
448  * which can be implemented as
449  *
450  * \f$ f_{pj} = \sum_{i=0}^{nq_0} \phi^A_p(\eta_{0,i})
451  * \tilde{u}_{i,j}
452  * \rightarrow {\bf B_1 U} \f$
453  * \f$ I_{pq} = \sum_{j=0}^{nq_1} \phi^B_{pq}(\eta_{1,j}) f_{pj}
454  * \rightarrow {\bf B_2[p*skip] f[skip]} \f$
455  *
456  * \b Recall: \f$ \eta_{1} = \frac{2(1+\xi_1)}{(1-\xi_2)}-1, \,
457  * \eta_2 = \xi_2\f$
458  *
459  * \b Note: For the orthgonality of this expansion to be realised the
460  * 'q' ordering must run fastest in contrast to the Quad and Hex
461  * ordering where 'p' index runs fastest to be consistent with the
462  * quadrature ordering.
463  *
464  * In the triangular space the i (i.e. \f$\eta_1\f$ direction)
465  * ordering still runs fastest by convention.
466  */
468  const Array<OneD, const NekDouble>& inarray,
469  Array<OneD, NekDouble>& outarray)
470  {
471  StdTriExp::v_IProductWRTBase_SumFac(inarray,outarray);
472  }
473 
475  const Array<OneD, const NekDouble>& inarray,
476  Array<OneD, NekDouble>& outarray)
477  {
478  int nq = GetTotPoints();
479  StdMatrixKey iprodmatkey(eIProductWRTBase,DetShapeType(),*this);
480  DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
481 
482  Blas::Dgemv('N',m_ncoeffs,nq,1.0,iprodmat->GetPtr().get(),
483  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
484  }
485 
487  const Array<OneD, const NekDouble>& inarray,
488  Array<OneD, NekDouble>& outarray,
489  bool multiplybyweights)
490  {
491  int nquad0 = m_base[0]->GetNumPoints();
492  int nquad1 = m_base[1]->GetNumPoints();
493  int order0 = m_base[0]->GetNumModes();
494 
495  if(multiplybyweights)
496  {
497  Array<OneD,NekDouble> tmp(nquad0*nquad1+nquad1*order0);
498  Array<OneD,NekDouble> wsp(tmp+nquad0*nquad1);
499 
500  // multiply by integration constants
501  MultiplyByQuadratureMetric(inarray,tmp);
502  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
503  m_base[1]->GetBdata(),
504  tmp,outarray,wsp);
505  }
506  else
507  {
508  Array<OneD,NekDouble> wsp(nquad1*order0);
509  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
510  m_base[1]->GetBdata(),
511  inarray,outarray,wsp);
512  }
513  }
514 
516  const Array<OneD, const NekDouble>& base0,
517  const Array<OneD, const NekDouble>& base1,
518  const Array<OneD, const NekDouble>& inarray,
519  Array<OneD, NekDouble>& outarray,
521  bool doCheckCollDir0,
522  bool doCheckCollDir1)
523  {
524  boost::ignore_unused(doCheckCollDir0, doCheckCollDir1);
525 
526  int i;
527  int mode;
528  int nquad0 = m_base[0]->GetNumPoints();
529  int nquad1 = m_base[1]->GetNumPoints();
530  int nmodes0 = m_base[0]->GetNumModes();
531  int nmodes1 = m_base[1]->GetNumModes();
532 
533  ASSERTL1(wsp.num_elements() >= nquad1*nmodes0,
534  "Workspace size is not sufficient");
535 
536  Blas::Dgemm('T','N',nquad1,nmodes0,nquad0,1.0,inarray.get(),nquad0,
537  base0.get(),nquad0,0.0,wsp.get(),nquad1);
538 
539  // Inner product with respect to 'b' direction
540  for (mode=i=0; i < nmodes0; ++i)
541  {
542  Blas::Dgemv('T',nquad1,nmodes1-i,1.0, base1.get()+mode*nquad1,
543  nquad1,wsp.get()+i*nquad1,1, 0.0,
544  outarray.get() + mode,1);
545  mode += nmodes1 - i;
546  }
547 
548  // fix for modified basis by splitting top vertex mode
550  {
551  outarray[1] += Blas::Ddot(nquad1,base1.get()+nquad1,1,
552  wsp.get()+nquad1,1);
553  }
554  }
555 
557  const int dir,
558  const Array<OneD, const NekDouble>& inarray,
559  Array<OneD, NekDouble>& outarray)
560  {
561  StdTriExp::v_IProductWRTDerivBase_SumFac(dir,inarray,outarray);
562  }
563 
565  const int dir,
566  const Array<OneD, const NekDouble>& inarray,
567  Array<OneD, NekDouble>& outarray)
568  {
569  int nq = GetTotPoints();
571 
572  switch(dir)
573  {
574  case 0:
575  {
576  mtype = eIProductWRTDerivBase0;
577  break;
578  }
579  case 1:
580  {
581  mtype = eIProductWRTDerivBase1;
582  break;
583  }
584  default:
585  {
586  ASSERTL1(false,"input dir is out of range");
587  break;
588  }
589  }
590 
591  StdMatrixKey iprodmatkey(mtype,DetShapeType(),*this);
592  DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
593 
594  Blas::Dgemv('N',m_ncoeffs,nq,1.0,iprodmat->GetPtr().get(),
595  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
596  }
597 
599  const int dir,
600  const Array<OneD, const NekDouble>& inarray,
601  Array<OneD, NekDouble>& outarray)
602  {
603  int i;
604  int nquad0 = m_base[0]->GetNumPoints();
605  int nquad1 = m_base[1]->GetNumPoints();
606  int nqtot = nquad0*nquad1;
607  int nmodes0 = m_base[0]->GetNumModes();
608  int wspsize = max(max(nqtot,m_ncoeffs),nquad1*nmodes0);
609 
610  Array<OneD, NekDouble> gfac0(2*wspsize);
611  Array<OneD, NekDouble> tmp0 (gfac0+wspsize);
612 
613  const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
614 
615  // set up geometric factor: 2/(1-z1)
616  for(i = 0; i < nquad1; ++i)
617  {
618  gfac0[i] = 2.0/(1-z1[i]);
619  }
620 
621  for(i = 0; i < nquad1; ++i)
622  {
623  Vmath::Smul(nquad0,gfac0[i],&inarray[0]+i*nquad0,1,
624  &tmp0[0]+i*nquad0,1);
625  }
626 
627  MultiplyByQuadratureMetric(tmp0,tmp0);
628 
629  switch(dir)
630  {
631  case 0:
632  {
633  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
634  m_base[1]->GetBdata(),
635  tmp0,outarray,gfac0);
636  break;
637  }
638  case 1:
639  {
641  const Array<OneD, const NekDouble>& z0 = m_base[0]->GetZ();
642 
643  for (i = 0; i < nquad0; ++i)
644  {
645  gfac0[i] = 0.5*(1+z0[i]);
646  }
647 
648  for (i = 0; i < nquad1; ++i)
649  {
650  Vmath::Vmul(nquad0,&gfac0[0],1,&tmp0[0]+i*nquad0,1,
651  &tmp0[0]+i*nquad0,1);
652  }
653 
654  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
655  m_base[1]->GetBdata(),
656  tmp0,tmp3,gfac0);
657 
658  MultiplyByQuadratureMetric(inarray,tmp0);
659  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
660  m_base[1]->GetDbdata(),
661  tmp0,outarray,gfac0);
662  Vmath::Vadd(m_ncoeffs,&tmp3[0],1,&outarray[0],1,
663  &outarray[0],1);
664  break;
665  }
666  default:
667  {
668  ASSERTL1(false, "input dir is out of range");
669  break;
670  }
671  }
672  }
673 
674  //---------------------------------------
675  // Evaluation functions
676  //---------------------------------------
677 
680  {
681 
682  // set up local coordinate system
683  if (fabs(xi[1]-1.0) < NekConstants::kNekZeroTol)
684  {
685  eta[0] = -1.0;
686  eta[1] = 1.0;
687  }
688  else
689  {
690  eta[0] = 2*(1+xi[0])/(1-xi[1])-1.0;
691  eta[1] = xi[1];
692  }
693  }
694 
696  const int mode, Array<OneD, NekDouble> &outarray)
697  {
698  int i,m;
699  int nquad0 = m_base[0]->GetNumPoints();
700  int nquad1 = m_base[1]->GetNumPoints();
701  int order0 = m_base[0]->GetNumModes();
702  int order1 = m_base[1]->GetNumModes();
703  int mode0 = 0;
704  Array<OneD, const NekDouble> base0 = m_base[0]->GetBdata();
705  Array<OneD, const NekDouble> base1 = m_base[1]->GetBdata();
706 
707  ASSERTL2(mode <= m_ncoeffs,
708  "calling argument mode is larger than "
709  "total expansion order");
710 
711  m = order1;
712  for (i = 0; i < order0; ++i, m+=order1-i)
713  {
714  if (m > mode)
715  {
716  mode0 = i;
717  break;
718  }
719  }
720 
721  // deal with top vertex mode in modified basis
722  if (mode == 1 &&
724  {
725  Vmath::Fill(nquad0*nquad1 , 1.0, outarray, 1);
726  }
727  else
728  {
729  for (i = 0; i < nquad1; ++i)
730  {
731  Vmath::Vcopy(nquad0,(NekDouble *)(base0.get()+mode0*nquad0),
732  1,&outarray[0]+i*nquad0,1);
733  }
734  }
735 
736  for(i = 0; i < nquad0; ++i)
737  {
738  Vmath::Vmul(nquad1,(NekDouble *)(base1.get() + mode*nquad1),
739  1,&outarray[0]+i,nquad0,&outarray[0]+i,nquad0);
740  }
741  }
742 
743 
745  {
746  return 3;
747  }
748 
750  {
751  return 3;
752  }
753 
755  {
757  }
758 
760  {
762  "BasisType is not a boundary interior form");
764  "BasisType is not a boundary interior form");
765 
766  return 3 + (GetBasisNumModes(0)-2) + 2*(GetBasisNumModes(1)-2);
767  }
768 
770  {
772  "BasisType is not a boundary interior form");
774  "BasisType is not a boundary interior form");
775 
776  return GetBasisNumModes(0) + 2*GetBasisNumModes(1);
777  }
778 
779  int StdTriExp::v_GetEdgeNcoeffs(const int i) const
780  {
781  ASSERTL2(i >= 0 && i <= 2, "edge id is out of range");
782 
783  if (i == 0)
784  {
785  return GetBasisNumModes(0);
786  }
787  else
788  {
789  return GetBasisNumModes(1);
790  }
791  }
792 
793  int StdTriExp::v_GetEdgeNumPoints(const int i) const
794  {
795  ASSERTL2((i >= 0)&&(i <= 2),"edge id is out of range");
796 
797  if (i == 0)
798  {
799  return GetNumPoints(0);
800  }
801  else
802  {
803  return GetNumPoints(1);
804  }
805  }
806 
808  const std::vector<unsigned int> &nummodes,
809  int &modes_offset)
810  {
812  nummodes[modes_offset],
813  nummodes[modes_offset+1]);
814  modes_offset += 2;
815 
816  return nmodes;
817  }
818 
820  {
821  ASSERTL2(i >= 0 && i <= 2, "edge id is out of range");
822 
823  if (i == 0)
824  {
825  return GetBasisType(0);
826  }
827  else
828  {
829  return GetBasisType(1);
830  }
831  }
832 
833 
835  Array<OneD, NekDouble> &coords_1,
836  Array<OneD, NekDouble> &coords_2)
837  {
838  boost::ignore_unused(coords_2);
839 
840  Array<OneD, const NekDouble> z0 = m_base[0]->GetZ();
841  Array<OneD, const NekDouble> z1 = m_base[1]->GetZ();
842  int nq0 = GetNumPoints(0);
843  int nq1 = GetNumPoints(1);
844  int i,j;
845 
846  for(i = 0; i < nq1; ++i)
847  {
848  for(j = 0; j < nq0; ++j)
849  {
850  coords_0[i*nq0+j] = (1+z0[j])*(1-z1[i])/2.0 - 1.0;
851  }
852  Vmath::Fill(nq0,z1[i],&coords_1[0] + i*nq0,1);
853  }
854  }
855 
857  {
858  return m_base[0]->GetBasisType() == LibUtilities::eModified_A &&
859  m_base[1]->GetBasisType() == LibUtilities::eModified_B;
860  }
861 
863  {
864  ASSERTL2(edge >= 0 && edge <= 2, "edge id is out of range");
865 
866  return edge == 0 ? 0 : 1;
867  }
868 
870  const int i) const
871  {
872  ASSERTL2(i >= 0 && i <= 2, "edge id is out of range");
873 
874  if (i == 0)
875  {
876  return GetBasis(0)->GetBasisKey();
877  }
878  else
879  {
880  switch(m_base[1]->GetBasisType())
881  {
883  {
884  switch(m_base[1]->GetPointsType())
885  {
887  {
889  m_base[1]->GetBasisKey().GetPointsKey().
890  GetNumPoints()+1,
892  return LibUtilities::BasisKey(
894  m_base[1]->GetNumModes(),pkey);
895  break;
896  }
897 
898  default:
900  "unexpected points distribution");
901  break;
902  }
903  break;
904  }
905  default:
907  "Information not available to set edge key");
908  break;
909  }
910  }
912  }
913 
914 
915 
916  //--------------------------
917  // Mappings
918  //--------------------------
919 
921  const int eid,
922  const Orientation edgeOrient,
923  Array<OneD, unsigned int>& maparray,
924  Array<OneD, int>& signarray,
925  int P)
926  {
929  "Mapping not defined for this type of basis");
930 
931  int i;
932  int numModes=0;
933  int order0 = m_base[0]->GetNumModes();
934  int order1 = m_base[1]->GetNumModes();
935 
936  switch (eid)
937  {
938  case 0:
939  numModes = order0;
940  break;
941  case 1:
942  case 2:
943  numModes = order1;
944  break;
945  default:
946  ASSERTL0(false,"eid must be between 0 and 2");
947  }
948 
949  bool checkForZeroedModes = false;
950  if (P == -1)
951  {
952  P = numModes;
953  }
954  else if(P != numModes)
955  {
956  checkForZeroedModes = true;
957  }
958 
959 
960  if(maparray.num_elements() != P)
961  {
962  maparray = Array<OneD, unsigned int>(P);
963  }
964 
965  if(signarray.num_elements() != P)
966  {
967  signarray = Array<OneD, int>(P,1);
968  }
969  else
970  {
971  fill(signarray.get() , signarray.get()+P, 1);
972  }
973 
974  switch(eid)
975  {
976  case 0:
977  {
978  int cnt = 0;
979  for(i = 0; i < P; cnt+=order1-i, ++i)
980  {
981  maparray[i] = cnt;
982  }
983 
984  if(edgeOrient==eBackwards)
985  {
986  swap( maparray[0] , maparray[1] );
987 
988  for(i = 3; i < P; i+=2)
989  {
990  signarray[i] = -1;
991  }
992  }
993  break;
994  }
995  case 1:
996  {
997  maparray[0] = order1;
998  maparray[1] = 1;
999  for(i = 2; i < P; i++)
1000  {
1001  maparray[i] = order1-1+i;
1002  }
1003 
1004  if(edgeOrient==eBackwards)
1005  {
1006  swap( maparray[0] , maparray[1] );
1007 
1008  for(i = 3; i < P; i+=2)
1009  {
1010  signarray[i] = -1;
1011  }
1012  }
1013  break;
1014  }
1015  case 2:
1016  {
1017  for(i = 0; i < P; i++)
1018  {
1019  maparray[i] = i;
1020  }
1021 
1022  if(edgeOrient==eBackwards)
1023  {
1024  swap( maparray[0] , maparray[1] );
1025 
1026  for(i = 3; i < P; i+=2)
1027  {
1028  signarray[i] = -1;
1029  }
1030  }
1031  break;
1032  }
1033  default:
1034  ASSERTL0(false,"eid must be between 0 and 2");
1035  break;
1036  }
1037 
1038 
1039  if (checkForZeroedModes)
1040  {
1041  // Zero signmap and set maparray to zero if
1042  // elemental modes are not as large as face modes
1043  for (int j = numModes; j < P; j++)
1044  {
1045  signarray[j] = 0.0;
1046  maparray[j] = maparray[0];
1047  }
1048  }
1049  }
1050 
1051  int StdTriExp::v_GetVertexMap(const int localVertexId,bool useCoeffPacking)
1052  {
1053  ASSERTL0(
1054  GetEdgeBasisType(localVertexId) == LibUtilities::eModified_A ||
1055  GetEdgeBasisType(localVertexId) == LibUtilities::eModified_B,
1056  "Mapping not defined for this type of basis");
1057 
1058  int localDOF = 0;
1059  if(useCoeffPacking == true)
1060  {
1061  switch(localVertexId)
1062  {
1063  case 0:
1064  {
1065  localDOF = 0;
1066  break;
1067  }
1068  case 1:
1069  {
1070  localDOF = 1;
1071  break;
1072  }
1073  case 2:
1074  {
1075  localDOF = m_base[1]->GetNumModes();
1076  break;
1077  }
1078  default:
1079  {
1080  ASSERTL0(false,"eid must be between 0 and 2");
1081  break;
1082  }
1083  }
1084  }
1085  else // follow book format for vertex indexing.
1086  {
1087  switch(localVertexId)
1088  {
1089  case 0:
1090  {
1091  localDOF = 0;
1092  break;
1093  }
1094  case 1:
1095  {
1096  localDOF = m_base[1]->GetNumModes();
1097  break;
1098  }
1099  case 2:
1100  {
1101  localDOF = 1;
1102  break;
1103  }
1104  default:
1105  {
1106  ASSERTL0(false,"eid must be between 0 and 2");
1107  break;
1108  }
1109  }
1110  }
1111 
1112  return localDOF;
1113  }
1114 
1116  const int eid,
1117  const Orientation edgeOrient,
1118  Array<OneD, unsigned int>& maparray,
1119  Array<OneD, int>& signarray)
1120  {
1123  "Mapping not defined for this type of basis");
1124  int i;
1125  const int nummodes1 = m_base[1]->GetNumModes();
1126  const int nEdgeIntCoeffs = GetEdgeNcoeffs(eid)-2;
1127 
1128  if(maparray.num_elements() != nEdgeIntCoeffs)
1129  {
1130  maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
1131  }
1132 
1133  if(signarray.num_elements() != nEdgeIntCoeffs)
1134  {
1135  signarray = Array<OneD, int>(nEdgeIntCoeffs,1);
1136  }
1137  else
1138  {
1139  fill( signarray.get() , signarray.get()+nEdgeIntCoeffs, 1 );
1140  }
1141 
1142  switch(eid)
1143  {
1144  case 0:
1145  {
1146  int cnt = 2*nummodes1 - 1;
1147  for(i = 0; i < nEdgeIntCoeffs; cnt+=nummodes1-2-i, ++i)
1148  {
1149  maparray[i] = cnt;
1150  }
1151 
1152  if(edgeOrient==eBackwards)
1153  {
1154  for(i = 1; i < nEdgeIntCoeffs; i+=2)
1155  {
1156  signarray[i] = -1;
1157  }
1158  }
1159  break;
1160  }
1161  case 1:
1162  {
1163  for(i = 0; i < nEdgeIntCoeffs; i++)
1164  {
1165  maparray[i] = nummodes1+1+i;
1166  }
1167 
1168  if(edgeOrient==eBackwards)
1169  {
1170  for(i = 1; i < nEdgeIntCoeffs; i+=2)
1171  {
1172  signarray[i] = -1;
1173  }
1174  }
1175  break;
1176  }
1177  case 2:
1178  {
1179  for(i = 0; i < nEdgeIntCoeffs; i++)
1180  {
1181  maparray[i] = 2+i;
1182  }
1183 
1184  if(edgeOrient==eBackwards)
1185  {
1186  for(i = 1; i < nEdgeIntCoeffs; i+=2)
1187  {
1188  signarray[i] = -1;
1189  }
1190  }
1191  break;
1192  }
1193  default:
1194  {
1195  ASSERTL0(false,"eid must be between 0 and 2");
1196  break;
1197  }
1198  }
1199  }
1200 
1202  {
1205  "Expansion not of a proper type");
1206 
1207  int i,j;
1208  int cnt = 0;
1209  int nummodes0, nummodes1;
1210  int startvalue;
1211  if(outarray.num_elements()!=GetNcoeffs()-NumBndryCoeffs())
1212  {
1214  }
1215 
1216  nummodes0 = m_base[0]->GetNumModes();
1217  nummodes1 = m_base[1]->GetNumModes();
1218 
1219  startvalue = 2*nummodes1;
1220 
1221  for(i = 0; i < nummodes0-2; i++)
1222  {
1223  for(j = 0; j < nummodes1-3-i; j++)
1224  {
1225  outarray[cnt++]=startvalue+j;
1226  }
1227  startvalue+=nummodes1-2-i;
1228  }
1229  }
1230 
1232  {
1235  "Expansion not of expected type");
1236  int i;
1237  int cnt;
1238  int nummodes0, nummodes1;
1239  int value;
1240 
1241  if (outarray.num_elements()!=NumBndryCoeffs())
1242  {
1244  }
1245 
1246  nummodes0 = m_base[0]->GetNumModes();
1247  nummodes1 = m_base[1]->GetNumModes();
1248 
1249  value = 2*nummodes1-1;
1250  for(i = 0; i < value; i++)
1251  {
1252  outarray[i]=i;
1253  }
1254  cnt = value;
1255 
1256  for(i = 0; i < nummodes0-2; i++)
1257  {
1258  outarray[cnt++]=value;
1259  value += nummodes1-2-i;
1260  }
1261  }
1262 
1263 
1264  //---------------------------------------
1265  // Wrapper functions
1266  //---------------------------------------
1267 
1269  {
1270 
1271  MatrixType mtype = mkey.GetMatrixType();
1272 
1273  DNekMatSharedPtr Mat;
1274 
1275  switch(mtype)
1276  {
1278  {
1279  int nq0, nq1, nq;
1280 
1281  nq0 = m_base[0]->GetNumPoints();
1282  nq1 = m_base[1]->GetNumPoints();
1283 
1284  // take definition from key
1286  {
1287  nq = (int) mkey.GetConstFactor(eFactorConst);
1288  }
1289  else
1290  {
1291  nq = max(nq0,nq1);
1292  }
1293 
1294  int neq = LibUtilities::StdTriData::
1296  Array<OneD, Array<OneD, NekDouble> > coords(neq);
1297  Array<OneD, NekDouble> coll (2);
1299  Array<OneD, NekDouble> tmp (nq0);
1300 
1301 
1302  Mat = MemoryManager<DNekMat>::AllocateSharedPtr(neq,nq0*nq1);
1303  int cnt = 0;
1304 
1305  for(int i = 0; i < nq; ++i)
1306  {
1307  for(int j = 0; j < nq-i; ++j,++cnt)
1308  {
1309  coords[cnt] = Array<OneD, NekDouble>(2);
1310  coords[cnt][0] = -1.0 + 2*j/(NekDouble)(nq-1);
1311  coords[cnt][1] = -1.0 + 2*i/(NekDouble)(nq-1);
1312  }
1313  }
1314 
1315  for(int i = 0; i < neq; ++i)
1316  {
1317  LocCoordToLocCollapsed(coords[i],coll);
1318 
1319  I[0] = m_base[0]->GetI(coll);
1320  I[1] = m_base[1]->GetI(coll+1);
1321 
1322  // interpolate first coordinate direction
1323  for (int j = 0; j < nq1; ++j)
1324  {
1325  NekDouble fac = (I[1]->GetPtr())[j];
1326  Vmath::Smul(nq0, fac, I[0]->GetPtr(), 1, tmp, 1);
1327 
1328  Vmath::Vcopy(nq0, &tmp[0], 1,
1329  Mat->GetRawPtr() + j*nq0*neq + i, neq);
1330  }
1331 
1332  }
1333  break;
1334  }
1335  default:
1336  {
1338  break;
1339  }
1340  }
1341 
1342  return Mat;
1343  }
1344 
1346  {
1347  return v_GenMatrix(mkey);
1348  }
1349 
1350 
1351  //---------------------------------------
1352  // Operator evaluation functions
1353  //---------------------------------------
1354 
1356  const Array<OneD, const NekDouble> &inarray,
1357  Array<OneD, NekDouble> &outarray,
1358  const StdMatrixKey &mkey)
1359  {
1360  StdExpansion::MassMatrixOp_MatFree(inarray,outarray,mkey);
1361  }
1362 
1364  const Array<OneD, const NekDouble> &inarray,
1365  Array<OneD, NekDouble> &outarray,
1366  const StdMatrixKey &mkey)
1367  {
1368  StdTriExp::v_LaplacianMatrixOp_MatFree(inarray,outarray,mkey);
1369  }
1370 
1372  const int k1,
1373  const int k2,
1374  const Array<OneD, const NekDouble> &inarray,
1375  Array<OneD, NekDouble> &outarray,
1376  const StdMatrixKey &mkey)
1377  {
1379  k1,k2,inarray,outarray,mkey);
1380  }
1381 
1383  const int i,
1384  const Array<OneD, const NekDouble> &inarray,
1385  Array<OneD, NekDouble> &outarray,
1386  const StdMatrixKey &mkey)
1387  {
1388  StdExpansion::WeakDerivMatrixOp_MatFree(i,inarray,outarray,mkey);
1389  }
1390 
1392  const Array<OneD, const NekDouble> &inarray,
1393  Array<OneD, NekDouble> &outarray,
1394  const StdMatrixKey &mkey)
1395  {
1396  StdTriExp::v_HelmholtzMatrixOp_MatFree(inarray,outarray,mkey);
1397  }
1398 
1399 
1401  const StdMatrixKey &mkey)
1402  {
1403  int qa = m_base[0]->GetNumPoints();
1404  int qb = m_base[1]->GetNumPoints();
1405  int nmodes_a = m_base[0]->GetNumModes();
1406  int nmodes_b = m_base[1]->GetNumModes();
1407 
1408  // Declare orthogonal basis.
1411 
1414  StdTriExp OrthoExp(Ba,Bb);
1415 
1416  Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
1417 
1418  // project onto physical space.
1419  OrthoExp.FwdTrans(array,orthocoeffs);
1420 
1421  if(mkey.ConstFactorExists(eFactorSVVPowerKerDiffCoeff)) // Rodrigo's power kern
1422  {
1424  NekDouble SvvDiffCoeff =
1427 
1428  int cnt = 0;
1429  for(int j = 0; j < nmodes_a; ++j)
1430  {
1431  for(int k = 0; k < nmodes_b-j; ++k, ++cnt)
1432  {
1433  NekDouble fac = std::max(
1434  pow((1.0*j)/(nmodes_a-1),cutoff*nmodes_a),
1435  pow((1.0*k)/(nmodes_b-1),cutoff*nmodes_b));
1436 
1437  orthocoeffs[cnt] *= (SvvDiffCoeff *fac);
1438  }
1439  }
1440  }
1441  else if(mkey.ConstFactorExists(eFactorSVVDGKerDiffCoeff)) // Rodrigo/mansoor's DG kernel
1442  {
1443  NekDouble SvvDiffCoeff =
1446  int max_ab = max(nmodes_a-kSVVDGFiltermodesmin,
1447  nmodes_b-kSVVDGFiltermodesmin);
1448  max_ab = max(max_ab,0);
1449  max_ab = min(max_ab,kSVVDGFiltermodesmax-kSVVDGFiltermodesmin);
1450 
1451  int cnt = 0;
1452  for(int j = 0; j < nmodes_a; ++j)
1453  {
1454  for(int k = 0; k < nmodes_b-j; ++k, ++cnt)
1455  {
1456  int maxjk = max(j,k);
1457  maxjk = min(maxjk,kSVVDGFiltermodesmax-1);
1458 
1459  orthocoeffs[cnt] *= SvvDiffCoeff *
1460  kSVVDGFilter[max_ab][maxjk];
1461  }
1462  }
1463  }
1464  else
1465  {
1466  NekDouble SvvDiffCoeff =
1468 
1469  int cutoff = (int) (mkey.GetConstFactor(eFactorSVVCutoffRatio)*
1470  min(nmodes_a,nmodes_b));
1471 
1472  NekDouble epsilon = 1.0;
1473  int nmodes = min(nmodes_a,nmodes_b);
1474 
1475  int cnt = 0;
1476 
1477  // apply SVV filter (JEL)
1478  for(int j = 0; j < nmodes_a; ++j)
1479  {
1480  for(int k = 0; k < nmodes_b-j; ++k)
1481  {
1482  if(j + k >= cutoff)
1483  {
1484  orthocoeffs[cnt] *= (SvvDiffCoeff
1485  *exp(-(j+k-nmodes)*(j+k-nmodes)
1486  /((NekDouble)((j+k-cutoff+epsilon)
1487  *(j+k-cutoff+epsilon)))));
1488  }
1489  else
1490  {
1491  orthocoeffs[cnt] *= 0.0;
1492  }
1493  cnt++;
1494  }
1495  }
1496 
1497  }
1498 
1499  // backward transform to physical space
1500  OrthoExp.BwdTrans(orthocoeffs,array);
1501  }
1502 
1504  int numMin,
1505  const Array<OneD, const NekDouble> &inarray,
1506  Array<OneD, NekDouble> &outarray)
1507  {
1508  int n_coeffs = inarray.num_elements();
1509  int nquad0 = m_base[0]->GetNumPoints();
1510  int nquad1 = m_base[1]->GetNumPoints();
1511  Array<OneD, NekDouble> coeff(n_coeffs);
1512  Array<OneD, NekDouble> coeff_tmp(n_coeffs,0.0);
1515  int nqtot = nquad0*nquad1;
1516  Array<OneD, NekDouble> phys_tmp(nqtot,0.0);
1517 
1518  int nmodes0 = m_base[0]->GetNumModes();
1519  int nmodes1 = m_base[1]->GetNumModes();
1520  int numMin2 = nmodes0;
1521  int i;
1522 
1523  const LibUtilities::PointsKey Pkey0(
1525  const LibUtilities::PointsKey Pkey1(
1527 
1528  LibUtilities::BasisKey b0(m_base[0]->GetBasisType(),nmodes0,Pkey0);
1529  LibUtilities::BasisKey b1(m_base[1]->GetBasisType(),nmodes1,Pkey1);
1530 
1531  LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A,nmodes0,Pkey0);
1532  LibUtilities::BasisKey bortho1(LibUtilities::eOrtho_B,nmodes1,Pkey1);
1533 
1534  StdRegions::StdTriExpSharedPtr m_OrthoTriExp;
1536 
1538  ::AllocateSharedPtr(b0, b1);
1539  m_OrthoTriExp = MemoryManager<StdRegions::StdTriExp>
1540  ::AllocateSharedPtr(bortho0, bortho1);
1541 
1542  m_TriExp ->BwdTrans(inarray,phys_tmp);
1543  m_OrthoTriExp->FwdTrans(phys_tmp, coeff);
1544 
1545  for (i = 0; i < n_coeffs; i++)
1546  {
1547  if (i == numMin)
1548  {
1549  coeff[i] = 0.0;
1550  numMin += numMin2 - 1;
1551  numMin2 -= 1.0;
1552  }
1553  }
1554 
1555  m_OrthoTriExp->BwdTrans(coeff,phys_tmp);
1556  m_TriExp ->FwdTrans(phys_tmp, outarray);
1557 
1558  }
1559 
1561  const Array<OneD, const NekDouble> &inarray,
1562  Array<OneD, NekDouble> &outarray,
1563  const StdMatrixKey &mkey)
1564  {
1566 
1567  if(inarray.get() == outarray.get())
1568  {
1570  Vmath::Vcopy(m_ncoeffs,inarray.get(),1,tmp.get(),1);
1571 
1572  Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, 1.0, mat->GetPtr().get(),
1573  m_ncoeffs, tmp.get(), 1, 0.0, outarray.get(), 1);
1574  }
1575  else
1576  {
1577  Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, 1.0, mat->GetPtr().get(),
1578  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
1579  }
1580  }
1581 
1582  //---------------------------------------
1583  // Private helper functions
1584  //---------------------------------------
1585 
1587  const Array<OneD, const NekDouble>& inarray,
1588  Array<OneD, NekDouble> &outarray)
1589  {
1590  int i;
1591  int nquad0 = m_base[0]->GetNumPoints();
1592  int nquad1 = m_base[1]->GetNumPoints();
1593 
1594  const Array<OneD, const NekDouble>& w0 = m_base[0]->GetW();
1595  const Array<OneD, const NekDouble>& w1 = m_base[1]->GetW();
1596  const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
1597 
1598  // multiply by integration constants
1599  for(i = 0; i < nquad1; ++i)
1600  {
1601  Vmath::Vmul(nquad0,inarray.get()+i*nquad0,1,
1602  w0.get(),1, outarray.get()+i*nquad0,1);
1603  }
1604 
1605  switch(m_base[1]->GetPointsType())
1606  {
1607  // Legendre inner product
1610  for(i = 0; i < nquad1; ++i)
1611  {
1612  Blas::Dscal(nquad0,0.5*(1-z1[i])*w1[i],
1613  outarray.get()+i*nquad0,1);
1614  }
1615  break;
1616 
1617  // (1,0) Jacobi Inner product
1619  for(i = 0; i < nquad1; ++i)
1620  {
1621  Blas::Dscal(nquad0,0.5*w1[i],outarray.get()+i*nquad0,1);
1622  }
1623  break;
1624 
1625  default:
1626  ASSERTL0(false, "Unsupported quadrature points type.");
1627  break;
1628  }
1629  }
1630 
1632  Array<OneD, int> &conn,
1633  bool standard)
1634  {
1635  boost::ignore_unused(standard);
1636 
1637  int np1 = m_base[0]->GetNumPoints();
1638  int np2 = m_base[1]->GetNumPoints();
1639  int np = max(np1,np2);
1640 
1641  conn = Array<OneD, int>(3*(np-1)*(np-1));
1642 
1643  int row = 0;
1644  int rowp1 = 0;
1645  int cnt = 0;
1646  for(int i = 0; i < np-1; ++i)
1647  {
1648  rowp1 += np-i;
1649  for(int j = 0; j < np-i-2; ++j)
1650  {
1651  conn[cnt++] = row +j;
1652  conn[cnt++] = row +j+1;
1653  conn[cnt++] = rowp1 +j;
1654 
1655  conn[cnt++] = rowp1 +j+1;
1656  conn[cnt++] = rowp1 +j;
1657  conn[cnt++] = row +j+1;
1658  }
1659 
1660  conn[cnt++] = row +np-i-2;
1661  conn[cnt++] = row +np-i-1;
1662  conn[cnt++] = rowp1+np-i-2;
1663 
1664  row += np-i;
1665  }
1666  }
1667 
1668 
1669  }//end namespace
1670 }//end namespace
virtual void v_IProductWRTDerivBase_MatOp(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdTriExp.cpp:564
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
virtual int v_GetEdgeNcoeffs(const int i) const
Definition: StdTriExp.cpp:779
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:209
virtual void v_GeneralMatrixOp_MatOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdTriExp.cpp:1560
const int kSVVDGFiltermodesmax
Definition: StdRegions.hpp:385
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Backward tranform for triangular elements.
Definition: StdTriExp.cpp:251
virtual void v_GetEdgeToElementMap(const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int P=-1)
Definition: StdTriExp.cpp:920
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
std::shared_ptr< StdTriExp > StdTriExpSharedPtr
Definition: StdTriExp.h:266
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:16
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
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:177
virtual void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdTriExp.cpp:1586
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_WeakDerivMatrixOp(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdTriExp.cpp:1382
void LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_FwdTrans_BndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdTriExp.cpp:331
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[p]*base1[pq] and put into ou...
Definition: StdTriExp.cpp:467
virtual void v_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdTriExp.cpp:1391
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:445
Principle Modified Functions .
Definition: BasisType.h:48
STL namespace.
void Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/y.
Definition: Vmath.cpp:274
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
NekDouble Integral(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &w0, const Array< OneD, const NekDouble > &w1)
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdTriExp.cpp:598
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
const int kSVVDGFiltermodesmin
Definition: StdRegions.hpp:384
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
LibUtilities::NekManager< StdMatrixKey, DNekBlkMat, StdMatrixKey::opLess > m_stdStaticCondMatrixManager
virtual int v_GetNverts() const
Definition: StdTriExp.cpp:744
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:81
virtual int v_NumBndryCoeffs() const
Definition: StdTriExp.cpp:759
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_GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdTriExp.cpp:1201
virtual void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
Definition: StdTriExp.cpp:1400
virtual const LibUtilities::BasisKey v_DetEdgeBasisKey(const int edge) const
Definition: StdTriExp.cpp:869
virtual void v_LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Definition: StdTriExp.cpp:678
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...
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:58
virtual int v_DetCartesianDirOfEdge(const int edge)
Definition: StdTriExp.cpp:862
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey)
Definition: StdTriExp.cpp:1268
Principle Orthogonal Functions .
Definition: BasisType.h:46
1D Evenly-spaced points using Lagrange polynomial
Definition: PointsType.h:64
virtual int v_CalcNumberOfCoefficients(const std::vector< unsigned int > &nummodes, int &modes_offset)
Definition: StdTriExp.cpp:807
virtual void v_GetSimplexEquiSpacedConnectivity(Array< OneD, int > &conn, bool standard=true)
Definition: StdTriExp.cpp:1631
virtual void v_ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdTriExp.cpp:1503
virtual int v_NumDGBndryCoeffs() const
Definition: StdTriExp.cpp:769
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
virtual void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray)
Definition: StdTriExp.cpp:1231
static const NekDouble kNekZeroTol
int GetNumModes() const
Returns the order of the basis.
Definition: Basis.h:86
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
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_MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdTriExp.cpp:1355
Principle Modified Functions .
Definition: BasisType.h:49
virtual void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdTriExp.cpp:1363
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey)
Definition: StdTriExp.cpp:1345
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.
int getNumberOfCoefficients(int Na, int Nb)
Definition: ShapeType.hpp:113
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
Principle Orthogonal Functions .
Definition: BasisType.h:45
Defines a specification for a set of points.
Definition: Points.h:59
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
Definition: StdTriExp.cpp:486
virtual void v_IProductWRTBase_MatOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdTriExp.cpp:474
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
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition: Vmath.cpp:318
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdTriExp.cpp:259
NekDouble GetConstFactor(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:121
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:130
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: StdTriExp.cpp:271
const NekDouble kSVVDGFilter[9][11]
Definition: StdRegions.hpp:387
virtual LibUtilities::ShapeType v_DetShapeType() const
Definition: StdTriExp.cpp:754
virtual int v_GetEdgeNumPoints(const int i) const
Definition: StdTriExp.cpp:793
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_StdPhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Definition: StdTriExp.cpp:223
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition: Blas.hpp:125
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray)
Integrates the specified function over the domain.
Definition: StdTriExp.cpp:84
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
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:215
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: StdTriExp.cpp:313
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdTriExp.cpp:556
bool ConstFactorExists(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:130
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:274
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:164
void MassMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual bool v_IsBoundaryInteriorExpansion()
Definition: StdTriExp.cpp:856
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_GetVertexMap(int localVertexId, bool useCoeffPacking=false)
Definition: StdTriExp.cpp:1051
virtual void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray)
Definition: StdTriExp.cpp:695
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.
LibUtilities::NekManager< StdMatrixKey, DNekMat, StdMatrixKey::opLess > m_stdMatrixManager
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140
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: StdTriExp.cpp:132
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_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: StdTriExp.cpp:515
std::shared_ptr< StdSegExp > StdSegExpSharedPtr
Definition: StdSegExp.h:46
virtual void v_GetEdgeInteriorMap(const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
Definition: StdTriExp.cpp:1115
#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
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
Definition: Blas.hpp:110
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
LibUtilities::BasisType GetEdgeBasisType(const int i) const
This function returns the type of expansion basis on the i-th edge.
Definition: StdExpansion.h:412
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_x, Array< OneD, NekDouble > &coords_y, Array< OneD, NekDouble > &coords_z)
Definition: StdTriExp.cpp:834
Describes the specification for a Basis.
Definition: Basis.h:49
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:302
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:186
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
virtual LibUtilities::BasisType v_GetEdgeBasisType(const int i) const
Definition: StdTriExp.cpp:819
virtual int v_GetNedges() const
Definition: StdTriExp.cpp:749
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Forward transformation from physical space to coefficient space...