Nektar++
TriExp.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File TriExp.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 // The above copyright notice and this permission notice shall be included
20 // in all copies or substantial portions of the Software.
21 //
22 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
23 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
24 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
25 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
26 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
27 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
28 // DEALINGS IN THE SOFTWARE.
29 //
30 // Description: Expasion for triangular elements.
31 //
32 ///////////////////////////////////////////////////////////////////////////////
33 
34 #include <boost/core/ignore_unused.hpp>
35 
37 #include <LocalRegions/TriExp.h>
38 #include <LocalRegions/SegExp.h>
42 
43 using namespace std;
44 
45 namespace Nektar
46 {
47  namespace LocalRegions
48  {
49  TriExp::TriExp(const LibUtilities::BasisKey &Ba,
50  const LibUtilities::BasisKey &Bb,
52  StdExpansion (LibUtilities::StdTriData::getNumberOfCoefficients(Ba.GetNumModes(),(Bb.GetNumModes())),2,Ba,Bb),
53  StdExpansion2D(LibUtilities::StdTriData::getNumberOfCoefficients(Ba.GetNumModes(),(Bb.GetNumModes())),Ba,Bb),
54  StdTriExp(Ba,Bb),
55  Expansion (geom),
56  Expansion2D (geom),
57  m_matrixManager(
58  std::bind(&TriExp::CreateMatrix, this, std::placeholders::_1),
59  std::string("TriExpMatrix")),
60  m_staticCondMatrixManager(
61  std::bind(&TriExp::CreateStaticCondMatrix, this, std::placeholders::_1),
62  std::string("TriExpStaticCondMatrix"))
63  {
64  }
65 
66 
68  StdExpansion(T),
69  StdExpansion2D(T),
70  StdTriExp(T),
71  Expansion(T),
72  Expansion2D(T),
73  m_matrixManager(T.m_matrixManager),
74  m_staticCondMatrixManager(T.m_staticCondMatrixManager)
75  {
76  }
77 
78 
80  {
81  }
82 
83 
84 
86  {
87  int nquad0 = m_base[0]->GetNumPoints();
88  int nquad1 = m_base[1]->GetNumPoints();
90  NekDouble ival;
91  Array<OneD,NekDouble> tmp(nquad0*nquad1);
92 
93  // multiply inarray with Jacobian
94  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
95  {
96  Vmath::Vmul(nquad0*nquad1, jac, 1, inarray, 1,tmp, 1);
97  }
98  else
99  {
100  Vmath::Smul(nquad0*nquad1, jac[0], inarray, 1, tmp, 1);
101  }
102 
103  // call StdQuadExp version;
104  ival = StdTriExp::v_Integral(tmp);
105  return ival;
106  }
107 
109  Array<OneD,NekDouble> &out_d0,
110  Array<OneD,NekDouble> &out_d1,
111  Array<OneD,NekDouble> &out_d2)
112  {
113  int nquad0 = m_base[0]->GetNumPoints();
114  int nquad1 = m_base[1]->GetNumPoints();
115  int nqtot = nquad0*nquad1;
117  = m_metricinfo->GetDerivFactors(GetPointsKeys());
118 
119  Array<OneD,NekDouble> diff0(2*nqtot);
120  Array<OneD,NekDouble> diff1(diff0+nqtot);
121 
122  StdTriExp::v_PhysDeriv(inarray, diff0, diff1);
123 
124  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
125  {
126  if(out_d0.size())
127  {
128  Vmath::Vmul (nqtot,df[0],1,diff0,1, out_d0, 1);
129  Vmath::Vvtvp (nqtot,df[1],1,diff1,1, out_d0, 1, out_d0,1);
130  }
131 
132  if(out_d1.size())
133  {
134  Vmath::Vmul (nqtot,df[2],1,diff0,1, out_d1, 1);
135  Vmath::Vvtvp (nqtot,df[3],1,diff1,1, out_d1, 1, out_d1,1);
136  }
137 
138  if(out_d2.size())
139  {
140  Vmath::Vmul (nqtot,df[4],1,diff0,1, out_d2, 1);
141  Vmath::Vvtvp (nqtot,df[5],1,diff1,1, out_d2, 1, out_d2,1);
142  }
143  }
144  else // regular geometry
145  {
146  if(out_d0.size())
147  {
148  Vmath::Smul (nqtot, df[0][0], diff0, 1, out_d0, 1);
149  Blas::Daxpy (nqtot, df[1][0], diff1, 1, out_d0, 1);
150  }
151 
152  if(out_d1.size())
153  {
154  Vmath::Smul (nqtot, df[2][0], diff0, 1, out_d1, 1);
155  Blas::Daxpy (nqtot, df[3][0], diff1, 1, out_d1, 1);
156  }
157 
158  if(out_d2.size())
159  {
160  Vmath::Smul (nqtot, df[4][0], diff0, 1, out_d2, 1);
161  Blas::Daxpy (nqtot, df[5][0], diff1, 1, out_d2, 1);
162  }
163  }
164  }
165 
166 
167  void TriExp::v_PhysDeriv(const int dir,
168  const Array<OneD, const NekDouble>& inarray,
169  Array<OneD, NekDouble> &outarray)
170  {
171  switch(dir)
172  {
173  case 0:
174  {
176  }
177  break;
178  case 1:
179  {
181  }
182  break;
183  case 2:
184  {
186  }
187  break;
188  default:
189  {
190  ASSERTL1(false,"input dir is out of range");
191  }
192  break;
193  }
194  }
195 
197  const Array<OneD, const NekDouble> &inarray,
198  const Array<OneD, const NekDouble> &direction,
200  {
201  if(! out.size())
202  {
203  return;
204  }
205 
206  int nquad0 = m_base[0]->GetNumPoints();
207  int nquad1 = m_base[1]->GetNumPoints();
208  int nqtot = nquad0*nquad1;
209 
210  const Array<TwoD, const NekDouble>& df =
211  m_metricinfo->GetDerivFactors(GetPointsKeys());
212 
213  Array<OneD,NekDouble> diff0(2*nqtot);
214  Array<OneD,NekDouble> diff1(diff0+nqtot);
215 
216  // diff0 = du/d_xi, diff1 = du/d_eta
217  StdTriExp::v_PhysDeriv(inarray, diff0, diff1);
218 
219  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
220  {
222 
223 
224  // D^v_xi = v_x*d_xi/dx + v_y*d_xi/dy + v_z*d_xi/dz
225  // D^v_eta = v_x*d_eta/dx + v_y*d_eta/dy + v_z*d_eta/dz
226  for (int i=0; i< 2; ++i)
227  {
228  tangmat[i] = Array<OneD, NekDouble>(nqtot,0.0);
229  for (int k=0; k<(m_geom->GetCoordim()); ++k)
230  {
231  Vmath::Vvtvp(nqtot,&df[2*k+i][0],1,&direction[k*nqtot],1,&tangmat[i][0],1,&tangmat[i][0],1);
232  }
233  }
234 
235  /// D_v = D^v_xi * du/d_xi + D^v_eta * du/d_eta
236  Vmath::Vmul (nqtot,&tangmat[0][0],1,&diff0[0],1, &out[0], 1);
237  Vmath::Vvtvp (nqtot,&tangmat[1][0],1,&diff1[0],1, &out[0], 1, &out[0],1);
238  }
239  else
240  {
242 
243  for (int i=0; i< 2; ++i)
244  {
245  tangmat[i] = Array<OneD, NekDouble>(nqtot,0.0);
246  for (int k=0; k<(m_geom->GetCoordim()); ++k)
247  {
248  Vmath::Svtvp(nqtot, df[2*k+i][0],
249  &direction[k*nqtot], 1,
250  &tangmat[i][0], 1, &tangmat[i][0], 1);
251  }
252  }
253 
254  /// D_v = D^v_xi * du/d_xi + D^v_eta * du/d_eta
255  Vmath::Vmul (nqtot, &tangmat[0][0], 1,
256  &diff0[0], 1,
257  &out[0], 1);
258 
259  Vmath::Vvtvp (nqtot, &tangmat[1][0], 1,
260  &diff1[0], 1,
261  &out[0], 1,
262  &out[0], 1);
263  }
264  }
265 
266 
268  Array<OneD,NekDouble> &outarray)
269  {
270  IProductWRTBase(inarray,outarray);
271 
272  // get Mass matrix inverse
274  DetShapeType(),*this);
275  DNekScalMatSharedPtr matsys = m_matrixManager[masskey];
276 
277  // copy inarray in case inarray == outarray
278  NekVector<NekDouble> in (m_ncoeffs,outarray,eCopy);
280 
281  out = (*matsys)*in;
282  }
283 
284 
286  Array<OneD, NekDouble> &outarray)
287  {
288  int i,j;
289  int npoints[2] = {m_base[0]->GetNumPoints(),
290  m_base[1]->GetNumPoints()};
291  int nmodes[2] = {m_base[0]->GetNumModes(),
292  m_base[1]->GetNumModes()};
293 
294  fill(outarray.get(), outarray.get()+m_ncoeffs, 0.0 );
295 
296  if(nmodes[0] == 1 && nmodes[1] == 1)
297  {
298  outarray[0] = inarray[0];
299  return;
300  }
301 
302  Array<OneD, NekDouble> physEdge[3];
303  Array<OneD, NekDouble> coeffEdge[3];
304  for(i = 0; i < 3; i++)
305  {
306  // define physEdge and add 1 so can interpolate grl10 points if necessary
307  physEdge[i] = Array<OneD, NekDouble>(max(npoints[i!=0],npoints[0]));
308  coeffEdge[i] = Array<OneD, NekDouble>(nmodes[i!=0]);
309  }
310 
311  for(i = 0; i < npoints[0]; i++)
312  {
313  physEdge[0][i] = inarray[i];
314  }
315 
316  // extract data in cartesian directions
317  for(i = 0; i < npoints[1]; i++)
318  {
319  physEdge[1][i] = inarray[npoints[0]-1+i*npoints[0]];
320  physEdge[2][i] = inarray[i*npoints[0]];
321  }
322 
323  SegExpSharedPtr segexp[3];
324  segexp[0] = MemoryManager<LocalRegions::SegExp>::AllocateSharedPtr(m_base[0]->GetBasisKey(),GetGeom2D()->GetEdge(0));
325 
327  {
328  for(i = 1; i < 3; i++)
329  {
330  segexp[i] = MemoryManager<LocalRegions::SegExp>::AllocateSharedPtr(m_base[i!=0]->GetBasisKey(),GetGeom2D()->GetEdge(i));
331  }
332  }
333  else // interploate using edge 0 GLL distribution
334  {
335  for(i = 1; i < 3; i++)
336  {
337  segexp[i] = MemoryManager<LocalRegions::SegExp>::AllocateSharedPtr(m_base[0]->GetBasisKey(),GetGeom2D()->GetEdge(i));
338 
339  LibUtilities::Interp1D(m_base[1]->GetPointsKey(),physEdge[i],
340  m_base[0]->GetPointsKey(),physEdge[i]);
341  }
342  npoints[1] = npoints[0];
343  }
344 
345 
346  Array<OneD, unsigned int> mapArray;
347  Array<OneD, int> signArray;
348  NekDouble sign;
349  // define an orientation to get EdgeToElmtMapping from Cartesian data
352 
353  for(i = 0; i < 3; i++)
354  {
355  segexp[i]->FwdTrans_BndConstrained(physEdge[i],coeffEdge[i]);
356 
357  // this orient goes with the one above and so could
358  // probably set both to eForwards
359  GetTraceToElementMap(i,mapArray,signArray,orient[i]);
360  for(j=0; j < nmodes[i!=0]; j++)
361  {
362  sign = (NekDouble) signArray[j];
363  outarray[ mapArray[j] ] = sign * coeffEdge[i][j];
364  }
365  }
366 
367  int nBoundaryDofs = NumBndryCoeffs();
368  int nInteriorDofs = m_ncoeffs - nBoundaryDofs;
369 
370  if (nInteriorDofs > 0) {
373 
375  MassMatrixOp(outarray,tmp0,stdmasskey);
376  IProductWRTBase(inarray,tmp1);
377 
378  Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
379 
380  // get Mass matrix inverse (only of interior DOF)
381  // use block (1,1) of the static condensed system
382  // note: this block alreay contains the inverse matrix
383  MatrixKey masskey(StdRegions::eMass,DetShapeType(),*this);
384  DNekScalMatSharedPtr matsys = (m_staticCondMatrixManager[masskey])->GetBlock(1,1);
385 
386  Array<OneD, NekDouble> rhs(nInteriorDofs);
387  Array<OneD, NekDouble> result(nInteriorDofs);
388 
389  GetInteriorMap(mapArray);
390 
391  for(i = 0; i < nInteriorDofs; i++)
392  {
393  rhs[i] = tmp1[ mapArray[i] ];
394  }
395 
396  Blas::Dgemv('N', nInteriorDofs, nInteriorDofs, matsys->Scale(), &((matsys->GetOwnedMatrix())->GetPtr())[0],
397  nInteriorDofs,rhs.get(),1,0.0,result.get(),1);
398 
399  for(i = 0; i < nInteriorDofs; i++)
400  {
401  outarray[ mapArray[i] ] = result[i];
402  }
403  }
404  }
405 
406 
408  Array<OneD, NekDouble> &outarray)
409  {
410  IProductWRTBase_SumFac(inarray,outarray);
411  }
412 
413 
414  void TriExp::v_IProductWRTDerivBase(const int dir,
415  const Array<OneD, const NekDouble>& inarray,
416  Array<OneD, NekDouble> & outarray)
417  {
418  IProductWRTDerivBase_SumFac(dir,inarray,outarray);
419  }
420 
421 
423  Array<OneD, NekDouble> &outarray,
424  bool multiplybyweights)
425  {
426  int nquad0 = m_base[0]->GetNumPoints();
427  int nquad1 = m_base[1]->GetNumPoints();
428  int order0 = m_base[0]->GetNumModes();
429 
430  if(multiplybyweights)
431  {
432  Array<OneD,NekDouble> tmp(nquad0*nquad1+nquad1*order0);
433  Array<OneD,NekDouble> wsp(tmp+nquad0*nquad1);
434 
435  MultiplyByQuadratureMetric(inarray,tmp);
436  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),m_base[1]->GetBdata(),tmp,outarray,wsp);
437  }
438  else
439  {
440  Array<OneD,NekDouble> wsp(+nquad1*order0);
441 
442  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),m_base[1]->GetBdata(),
443  inarray,outarray,wsp);
444  }
445  }
446 
447 
449  Array<OneD, NekDouble> &outarray)
450  {
451  int nq = GetTotPoints();
453  DNekScalMatSharedPtr iprodmat = m_matrixManager[iprodmatkey];
454 
455  Blas::Dgemv('N',m_ncoeffs,nq,iprodmat->Scale(),(iprodmat->GetOwnedMatrix())->GetPtr().get(),
456  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
457 
458  }
459 
460 
462  const Array<OneD, const NekDouble>& inarray,
463  Array<OneD, NekDouble> & outarray)
464  {
465  int nquad0 = m_base[0]->GetNumPoints();
466  int nquad1 = m_base[1]->GetNumPoints();
467  int nqtot = nquad0*nquad1;
468  int nmodes0 = m_base[0]->GetNumModes();
469  int wspsize = max(max(nqtot,m_ncoeffs),nquad1*nmodes0);
470 
471  Array<OneD, NekDouble> tmp0 (4*wspsize);
472  Array<OneD, NekDouble> tmp1 (tmp0 + wspsize);
473  Array<OneD, NekDouble> tmp2 (tmp0 + 2*wspsize);
474  Array<OneD, NekDouble> tmp3 (tmp0 + 3*wspsize);
475 
477  tmp2D[0] = tmp1;
478  tmp2D[1] = tmp2;
479 
480  TriExp::v_AlignVectorToCollapsedDir(dir, inarray, tmp2D);
481 
482  MultiplyByQuadratureMetric(tmp1,tmp1);
483  MultiplyByQuadratureMetric(tmp2,tmp2);
484 
485  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
486  m_base[1]->GetBdata() ,
487  tmp1,tmp3,tmp0);
488  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata() ,
489  m_base[1]->GetDbdata(),
490  tmp2,outarray,tmp0);
491  Vmath::Vadd(m_ncoeffs, tmp3, 1, outarray, 1, outarray, 1);
492  }
493 
495  const int dir,
496  const Array<OneD, const NekDouble> &inarray,
497  Array<OneD, Array<OneD, NekDouble> > &outarray)
498  {
499  ASSERTL1((dir==0)||(dir==1)||(dir==2),"Invalid direction.");
500  ASSERTL1((dir==2)?(m_geom->GetCoordim()==3):true,
501  "Invalid direction.");
502 
503  int nquad0 = m_base[0]->GetNumPoints();
504  int nquad1 = m_base[1]->GetNumPoints();
505  int nqtot = nquad0*nquad1;
506  int nmodes0 = m_base[0]->GetNumModes();
507  int wspsize = max(max(nqtot,m_ncoeffs),nquad1*nmodes0);
508 
509  const Array<TwoD, const NekDouble>& df =
510  m_metricinfo->GetDerivFactors(GetPointsKeys());
511 
512  Array<OneD, NekDouble> tmp0 (wspsize);
513  Array<OneD, NekDouble> tmp3 (wspsize);
514  Array<OneD, NekDouble> gfac0(wspsize);
515  Array<OneD, NekDouble> gfac1(wspsize);
516 
517  Array<OneD, NekDouble> tmp1 = outarray[0];
518  Array<OneD, NekDouble> tmp2 = outarray[1];
519 
520  const Array<OneD, const NekDouble>& z0 = m_base[0]->GetZ();
521  const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
522 
523  // set up geometric factor: 2/(1-z1)
524  for (int i = 0; i < nquad1; ++i)
525  {
526  gfac0[i] = 2.0/(1-z1[i]);
527  }
528  for (int i = 0; i < nquad0; ++i)
529  {
530  gfac1[i] = 0.5*(1+z0[i]);
531  }
532 
533  for (int i = 0; i < nquad1; ++i)
534  {
535  Vmath::Smul(nquad0, gfac0[i], &inarray[0]+i*nquad0, 1,
536  &tmp0[0]+i*nquad0, 1);
537  }
538 
539  for (int i = 0; i < nquad1; ++i)
540  {
541  Vmath::Vmul(nquad0, &gfac1[0], 1, &tmp0[0]+i*nquad0, 1,
542  &tmp1[0]+i*nquad0, 1);
543  }
544 
545  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
546  {
547  Vmath::Vmul(nqtot,&df[2*dir][0], 1,&tmp0[0], 1,&tmp0[0],1);
548  Vmath::Vmul(nqtot,&df[2*dir+1][0],1,&tmp1[0], 1,&tmp1[0],1);
549  Vmath::Vmul(nqtot,&df[2*dir+1][0],1,&inarray[0],1,&tmp2[0],1);
550  }
551  else
552  {
553  Vmath::Smul(nqtot, df[2*dir][0], tmp0, 1, tmp0, 1);
554  Vmath::Smul(nqtot, df[2*dir+1][0], tmp1, 1, tmp1, 1);
555  Vmath::Smul(nqtot, df[2*dir+1][0], inarray, 1, tmp2, 1);
556  }
557  Vmath::Vadd(nqtot, tmp0, 1, tmp1, 1, tmp1, 1);
558  }
559 
561  const Array<OneD, const NekDouble>& inarray,
562  Array<OneD, NekDouble> &outarray)
563  {
564  int nq = GetTotPoints();
566 
567  switch(dir)
568  {
569  case 0:
570  {
572  }
573  break;
574  case 1:
575  {
577  }
578  break;
579  case 2:
580  {
582  }
583  break;
584  default:
585  {
586  ASSERTL1(false,"input dir is out of range");
587  }
588  break;
589  }
590 
591  MatrixKey iprodmatkey(mtype,DetShapeType(),*this);
592  DNekScalMatSharedPtr iprodmat = m_matrixManager[iprodmatkey];
593 
594  Blas::Dgemv('N',m_ncoeffs,nq,iprodmat->Scale(),(iprodmat->GetOwnedMatrix())->GetPtr().get(),
595  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
596 
597  }
598 
599 
601  const Array<OneD, const NekDouble>& direction,
602  const Array<OneD, const NekDouble>& inarray,
603  Array<OneD, NekDouble> & outarray)
604  {
606  inarray,outarray);
607  }
608 
609 
610  /**
611  * @brief Directinoal Derivative in the modal space in the dir
612  * direction of varcoeffs.
613  */
615  const Array<OneD, const NekDouble>& direction,
616  const Array<OneD, const NekDouble>& inarray,
617  Array<OneD, NekDouble> & outarray)
618  {
619  int i;
620  int shapedim = 2;
621  int nquad0 = m_base[0]->GetNumPoints();
622  int nquad1 = m_base[1]->GetNumPoints();
623  int nqtot = nquad0*nquad1;
624  int nmodes0 = m_base[0]->GetNumModes();
625  int wspsize = max(max(nqtot,m_ncoeffs),nquad1*nmodes0);
626 
627  const Array<TwoD, const NekDouble>& df =
628  m_metricinfo->GetDerivFactors(GetPointsKeys());
629 
630  Array<OneD, NekDouble> tmp0 (6*wspsize);
631  Array<OneD, NekDouble> tmp1 (tmp0 + wspsize);
632  Array<OneD, NekDouble> tmp2 (tmp0 + 2*wspsize);
633  Array<OneD, NekDouble> tmp3 (tmp0 + 3*wspsize);
634  Array<OneD, NekDouble> gfac0(tmp0 + 4*wspsize);
635  Array<OneD, NekDouble> gfac1(tmp0 + 5*wspsize);
636 
637  const Array<OneD, const NekDouble>& z0 = m_base[0]->GetZ();
638  const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
639 
640  // set up geometric factor: 2/(1-z1)
641  for(i = 0; i < nquad1; ++i)
642  {
643  gfac0[i] = 2.0/(1-z1[i]);
644  }
645  for(i = 0; i < nquad0; ++i)
646  {
647  gfac1[i] = 0.5*(1+z0[i]);
648  }
649  for(i = 0; i < nquad1; ++i)
650  {
651  Vmath::Smul(nquad0, gfac0[i], &inarray[0] + i*nquad0, 1,
652  &tmp0[0] + i*nquad0, 1);
653  }
654  for(i = 0; i < nquad1; ++i)
655  {
656  Vmath::Vmul(nquad0, &gfac1[0], 1,
657  &tmp0[0] + i*nquad0, 1,
658  &tmp1[0] + i*nquad0, 1);
659  }
660 
661  // Compute gmat \cdot e^j
662  Array<OneD, Array<OneD, NekDouble> > dfdir(shapedim);
663  Expansion::ComputeGmatcdotMF(df, direction, dfdir);
664 
665  Vmath::Vmul(nqtot, &dfdir[0][0], 1, &tmp0[0], 1, &tmp0[0], 1);
666  Vmath::Vmul(nqtot, &dfdir[1][0], 1, &tmp1[0], 1, &tmp1[0], 1);
667  Vmath::Vmul(nqtot, &dfdir[1][0], 1, &inarray[0], 1, &tmp2[0], 1);
668 
669  Vmath::Vadd(nqtot, &tmp0[0], 1, &tmp1[0], 1, &tmp1[0], 1);
670 
671  MultiplyByQuadratureMetric(tmp1,tmp1);
672  MultiplyByQuadratureMetric(tmp2,tmp2);
673 
674  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
675  m_base[1]->GetBdata(),
676  tmp1, tmp3, tmp0);
677  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
678  m_base[1]->GetDbdata(),
679  tmp2, outarray, tmp0);
680  Vmath::Vadd(m_ncoeffs, tmp3, 1, outarray, 1, outarray, 1);
681  }
682 
683 
688  Array<OneD, NekDouble> &outarray)
689  {
690  int nq = m_base[0]->GetNumPoints()*m_base[1]->GetNumPoints();
692 
693  const Array<OneD, const Array<OneD, NekDouble> > &normals =
694  GetLeftAdjacentElementExp()->GetTraceNormal(
696 
697  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
698  {
699  Vmath::Vvtvvtp(nq,&normals[0][0],1,&Fx[0],1,
700  &normals[1][0],1,&Fy[0],1,&Fn[0],1);
701  Vmath::Vvtvp (nq,&normals[2][0],1,&Fz[0],1,&Fn[0],1,&Fn[0],1);
702  }
703  else
704  {
705  Vmath::Svtsvtp(nq,normals[0][0],&Fx[0],1,
706  normals[1][0],&Fy[0],1,&Fn[0],1);
707  Vmath::Svtvp (nq,normals[2][0],&Fz[0],1,&Fn[0],1,&Fn[0],1);
708  }
709 
710  IProductWRTBase(Fn,outarray);
711  }
712 
714  const Array<OneD, const Array<OneD, NekDouble> > &Fvec,
715  Array<OneD, NekDouble> &outarray)
716  {
717  NormVectorIProductWRTBase(Fvec[0], Fvec[1], Fvec[2], outarray);
718  }
719 
721  {
722 
724  ::AllocateSharedPtr(m_base[0]->GetBasisKey(),
725  m_base[1]->GetBasisKey());
726  }
727 
729  {
731  2, m_base[0]->GetPointsKey());
733  2, m_base[1]->GetPointsKey());
734 
736  ::AllocateSharedPtr( bkey0, bkey1);
737  }
738 
740  Array<OneD,NekDouble> &coords)
741  {
742  int i;
743 
744  ASSERTL1(Lcoords[0] >= -1.0 && Lcoords[1] <= 1.0 &&
745  Lcoords[1] >= -1.0 && Lcoords[1] <=1.0,
746  "Local coordinates are not in region [-1,1]");
747 
748  m_geom->FillGeom();
749 
750  for(i = 0; i < m_geom->GetCoordim(); ++i)
751  {
752  coords[i] = m_geom->GetCoord(i,Lcoords);
753  }
754  }
755 
757  Array<OneD, NekDouble> &coords_0,
758  Array<OneD, NekDouble> &coords_1,
759  Array<OneD, NekDouble> &coords_2)
760  {
761  Expansion::v_GetCoords(coords_0, coords_1, coords_2);
762  }
763 
764 
765  /**
766  * Given the local cartesian coordinate \a Lcoord evaluate the
767  * value of physvals at this point by calling through to the
768  * StdExpansion method
769  */
771  const Array<OneD, const NekDouble> &Lcoord,
772  const Array<OneD, const NekDouble> &physvals)
773  {
774  // Evaluate point in local (eta) coordinates.
775  return StdTriExp::v_PhysEvaluate(Lcoord,physvals);
776  }
777 
779  {
781 
782  ASSERTL0(m_geom,"m_geom not defined");
783  m_geom->GetLocCoords(coord,Lcoord);
784 
785  return StdTriExp::v_PhysEvaluate(Lcoord, physvals);
786  }
787 
788 
790  const int edge,
791  const StdRegions::StdExpansionSharedPtr &EdgeExp,
792  const Array<OneD, const NekDouble> &inarray,
793  Array<OneD,NekDouble> &outarray,
795  {
796  int nquad0 = m_base[0]->GetNumPoints();
797  int nquad1 = m_base[1]->GetNumPoints();
798 
799  // Extract in Cartesian direction because we have to deal with
800  // e.g. Gauss-Radau points.
801  switch(edge)
802  {
803  case 0:
804  Vmath::Vcopy(nquad0, &(inarray[0]), 1, &(outarray[0]), 1);
805  break;
806  case 1:
807  Vmath::Vcopy(nquad1, &(inarray[0])+(nquad0-1),
808  nquad0, &(outarray[0]), 1);
809  break;
810  case 2:
811  Vmath::Vcopy(nquad1, &(inarray[0]), nquad0, &(outarray[0]), 1);
812  break;
813  default:
814  ASSERTL0(false,"edge value (< 3) is out of range");
815  break;
816  }
817 
818  ASSERTL1(EdgeExp->GetBasis(0)->GetPointsType()
820  "Edge expansion should be GLL");
821 
822  // Interpolate if required
823  if(m_base[edge?1:0]->GetPointsKey() != EdgeExp->GetBasis(0)->GetPointsKey())
824  {
825  Array<OneD,NekDouble> outtmp(max(nquad0,nquad1));
826 
827  outtmp = outarray;
828 
829  LibUtilities::Interp1D(m_base[edge?1:0]->GetPointsKey(),
830  outtmp,
831  EdgeExp->GetBasis(0)->GetPointsKey(),
832  outarray);
833  }
834 
835  if (orient == StdRegions::eNoOrientation)
836  {
837  orient = GetTraceOrient(edge);
838  }
839 
840  //Reverse data if necessary
841  if(orient == StdRegions::eBackwards)
842  {
843  Vmath::Reverse(EdgeExp->GetNumPoints(0),&outarray[0],1,
844  &outarray[0],1);
845  }
846 
847  }
848 
849 
851  const int edge,const Array<OneD, const NekDouble> &inarray,
852  Array<OneD, NekDouble> &outarray)
853  {
854  boost::ignore_unused(edge, inarray, outarray);
855  ASSERTL0(false,
856  "Routine not implemented for triangular elements");
857  }
858 
860  const int edge,
861  Array<OneD, NekDouble> &outarray)
862  {
863  boost::ignore_unused(edge, outarray);
864  ASSERTL0(false,
865  "Routine not implemented for triangular elements");
866  }
867 
868 
869 
871  const int edge,
872  Array<OneD, int> &outarray)
873  {
874  int nquad0 = m_base[0]->GetNumPoints();
875  int nquad1 = m_base[1]->GetNumPoints();
876 
877  // Get points in Cartesian orientation
878  switch (edge)
879  {
880  case 0:
881  outarray = Array<OneD, int>(nquad0);
882  for (int i = 0; i < nquad0; ++i)
883  {
884  outarray[i] = i;
885  }
886  break;
887  case 1:
888  outarray = Array<OneD, int>(nquad1);
889  for (int i = 0; i < nquad1; ++i)
890  {
891  outarray[i] = (nquad0-1) + i * nquad0;
892  }
893  break;
894  case 2:
895  outarray = Array<OneD, int>(nquad1);
896  for (int i = 0; i < nquad1; ++i)
897  {
898  outarray[i] = i*nquad0;
899  }
900  break;
901  default:
902  ASSERTL0(false, "edge value (< 3) is out of range");
903  break;
904  }
905 
906  }
907 
908 
909  void TriExp::v_ComputeTraceNormal(const int edge)
910  {
911  int i;
912  const SpatialDomains::GeomFactorsSharedPtr & geomFactors = GetGeom()->GetMetricInfo();
913 
915  for(i = 0; i < ptsKeys.size(); ++i)
916  {
917  // Need at least 2 points for computing normals
918  if (ptsKeys[i].GetNumPoints() == 1)
919  {
920  LibUtilities::PointsKey pKey(2, ptsKeys[i].GetPointsType());
921  ptsKeys[i] = pKey;
922  }
923  }
924 
925  const SpatialDomains::GeomType type = geomFactors->GetGtype();
926  const Array<TwoD, const NekDouble> & df = geomFactors->GetDerivFactors(ptsKeys);
927  const Array<OneD, const NekDouble> & jac = geomFactors->GetJac(ptsKeys);
928  int nqe = m_base[0]->GetNumPoints();
929  int dim = GetCoordim();
930 
933  for (i = 0; i < dim; ++i)
934  {
935  normal[i] = Array<OneD, NekDouble>(nqe);
936  }
937 
938  size_t nqb = nqe;
939  size_t nbnd= edge;
942 
943  // Regular geometry case
945  {
946  NekDouble fac;
947  // Set up normals
948  switch(edge)
949  {
950  case 0:
951  for(i = 0; i < GetCoordim(); ++i)
952  {
953  Vmath::Fill(nqe,-df[2*i+1][0],normal[i],1);
954  }
955  break;
956  case 1:
957  for(i = 0; i < GetCoordim(); ++i)
958  {
959  Vmath::Fill(nqe,df[2*i+1][0] + df[2*i][0],normal[i],1);
960  }
961  break;
962  case 2:
963  for(i = 0; i < GetCoordim(); ++i)
964  {
965  Vmath::Fill(nqe,-df[2*i][0],normal[i],1);
966  }
967  break;
968  default:
969  ASSERTL0(false,"Edge is out of range (edge < 3)");
970  }
971 
972  // normalise
973  fac = 0.0;
974  for(i =0 ; i < GetCoordim(); ++i)
975  {
976  fac += normal[i][0]*normal[i][0];
977  }
978  fac = 1.0/sqrt(fac);
979 
980  Vmath::Fill(nqb, fac, length, 1);
981 
982  for (i = 0; i < GetCoordim(); ++i)
983  {
984  Vmath::Smul(nqe,fac,normal[i],1,normal[i],1);
985  }
986  }
987  else // Set up deformed normals
988  {
989  int j;
990 
991  int nquad0 = ptsKeys[0].GetNumPoints();
992  int nquad1 = ptsKeys[1].GetNumPoints();
993 
994  LibUtilities::PointsKey from_key;
995 
996  Array<OneD,NekDouble> normals(GetCoordim()*max(nquad0,nquad1),0.0);
997  Array<OneD,NekDouble> edgejac(GetCoordim()*max(nquad0,nquad1),0.0);
998 
999  // Extract Jacobian along edges and recover local
1000  // derivates (dx/dr) for polynomial interpolation by
1001  // multiplying m_gmat by jacobian
1002  switch(edge)
1003  {
1004  case 0:
1005  for(j = 0; j < nquad0; ++j)
1006  {
1007  edgejac[j] = jac[j];
1008  for(i = 0; i < GetCoordim(); ++i)
1009  {
1010  normals[i*nquad0+j] = -df[2*i+1][j]*edgejac[j];
1011  }
1012  }
1013  from_key = ptsKeys[0];
1014  break;
1015  case 1:
1016  for(j = 0; j < nquad1; ++j)
1017  {
1018  edgejac[j] = jac[nquad0*j+nquad0-1];
1019  for(i = 0; i < GetCoordim(); ++i)
1020  {
1021  normals[i*nquad1+j] = (df[2*i][nquad0*j + nquad0-1] + df[2*i+1][nquad0*j + nquad0-1])*edgejac[j];
1022  }
1023  }
1024  from_key = ptsKeys[1];
1025  break;
1026  case 2:
1027  for(j = 0; j < nquad1; ++j)
1028  {
1029  edgejac[j] = jac[nquad0*j];
1030  for(i = 0; i < GetCoordim(); ++i)
1031  {
1032  normals[i*nquad1+j] = -df[2*i][nquad0*j]*edgejac[j];
1033  }
1034  }
1035  from_key = ptsKeys[1];
1036  break;
1037  default:
1038  ASSERTL0(false,"edge is out of range (edge < 3)");
1039 
1040  }
1041 
1042  int nq = from_key.GetNumPoints();
1043  Array<OneD,NekDouble> work(nqe,0.0);
1044 
1045  // interpolate Jacobian and invert
1046  LibUtilities::Interp1D(from_key,jac,m_base[0]->GetPointsKey(),work);
1047  Vmath::Sdiv(nqe,1.0,&work[0],1,&work[0],1);
1048 
1049  // interpolate
1050  for(i = 0; i < GetCoordim(); ++i)
1051  {
1052  LibUtilities::Interp1D(from_key,&normals[i*nq],m_base[0]->GetPointsKey(),&normal[i][0]);
1053  Vmath::Vmul(nqe,work,1,normal[i],1,normal[i],1);
1054  }
1055 
1056  //normalise normal vectors
1057  Vmath::Zero(nqe,work,1);
1058  for(i = 0; i < GetCoordim(); ++i)
1059  {
1060  Vmath::Vvtvp(nqe,normal[i],1, normal[i],1,work,1,work,1);
1061  }
1062 
1063  Vmath::Vsqrt(nqe,work,1,work,1);
1064  Vmath::Sdiv(nqe,1.0,work,1,work,1);
1065 
1066  Vmath::Vcopy(nqb, work, 1, length, 1);
1067 
1068  for(i = 0; i < GetCoordim(); ++i)
1069  {
1070  Vmath::Vmul(nqe,normal[i],1,work,1,normal[i],1);
1071  }
1072  }
1073  if(GetGeom()->GetEorient(edge) == StdRegions::eBackwards)
1074  {
1075  for(i = 0; i < GetCoordim(); ++i)
1076  {
1077  if(geomFactors->GetGtype() == SpatialDomains::eDeformed)
1078  {
1079  Vmath::Reverse(nqe, normal[i], 1, normal[i],1);
1080  }
1081  }
1082  }
1083  }
1084 
1086  {
1087  return m_geom->GetCoordim();
1088  }
1089 
1090 
1092  const NekDouble *data,
1093  const std::vector<unsigned int > &nummodes,
1094  const int mode_offset,
1095  NekDouble * coeffs,
1096  std::vector<LibUtilities::BasisType> &fromType)
1097  {
1098  boost::ignore_unused(fromType);
1099 
1100  int data_order0 = nummodes[mode_offset];
1101  int fillorder0 = min(m_base[0]->GetNumModes(),data_order0);
1102  int data_order1 = nummodes[mode_offset+1];
1103  int order1 = m_base[1]->GetNumModes();
1104  int fillorder1 = min(order1,data_order1);
1105 
1106  switch(m_base[0]->GetBasisType())
1107  {
1110  {
1111  int i;
1112  int cnt = 0;
1113  int cnt1 = 0;
1114 
1117  "Extraction routine not set up for this basis");
1118 
1119  Vmath::Zero(m_ncoeffs,coeffs,1);
1120  for(i = 0; i < fillorder0; ++i)
1121  {
1122  Vmath::Vcopy(fillorder1-i,&data[cnt],1,&coeffs[cnt1],1);
1123  cnt += data_order1-i;
1124  cnt1 += order1-i;
1125  }
1126  }
1127  break;
1128  default:
1129  ASSERTL0(false,"basis is either not set up or not hierarchicial");
1130  }
1131  }
1132 
1133 
1135  {
1136  return GetGeom2D()->GetEorient(edge);
1137  }
1138 
1140  {
1141  ASSERTL1(dir >= 0 &&dir <= 1,"input dir is out of range");
1142  return m_base[dir];
1143  }
1144 
1145 
1146  int TriExp::v_GetNumPoints(const int dir) const
1147  {
1148  return GetNumPoints(dir);
1149  }
1150 
1151 
1153  {
1154  DNekMatSharedPtr returnval;
1155  switch(mkey.GetMatrixType())
1156  {
1164  returnval = Expansion2D::v_GenMatrix(mkey);
1165  break;
1166  default:
1167  returnval = StdTriExp::v_GenMatrix(mkey);
1168  break;
1169  }
1170 
1171  return returnval;
1172  }
1173 
1174 
1176  {
1177  LibUtilities::BasisKey bkey0 = m_base[0]->GetBasisKey();
1178  LibUtilities::BasisKey bkey1 = m_base[1]->GetBasisKey();
1180  AllocateSharedPtr(bkey0,bkey1);
1181 
1182  return tmp->GetStdMatrix(mkey);
1183  }
1184 
1185 
1187  {
1188  DNekScalMatSharedPtr returnval;
1190 
1191  ASSERTL2(m_metricinfo->GetGtype() != SpatialDomains::eNoGeomType,"Geometric information is not set up");
1192 
1193  switch(mkey.GetMatrixType())
1194  {
1195  case StdRegions::eMass:
1196  {
1197  if((m_metricinfo->GetGtype() == SpatialDomains::eDeformed)||
1198  (mkey.GetNVarCoeff()))
1199  {
1200  NekDouble one = 1.0;
1201  DNekMatSharedPtr mat = GenMatrix(mkey);
1202  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,mat);
1203  }
1204  else
1205  {
1206  NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
1207  DNekMatSharedPtr mat = GetStdMatrix(mkey);
1208  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(jac,mat);
1209  }
1210  }
1211  break;
1212  case StdRegions::eInvMass:
1213  {
1214  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1215  {
1216  NekDouble one = 1.0;
1218  *this);
1219  DNekMatSharedPtr mat = GenMatrix(masskey);
1220  mat->Invert();
1221 
1222  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,mat);
1223  }
1224  else
1225  {
1226  NekDouble fac = 1.0/(m_metricinfo->GetJac(ptsKeys))[0];
1227  DNekMatSharedPtr mat = GetStdMatrix(mkey);
1228  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(fac,mat);
1229 
1230  }
1231  }
1232  break;
1236  {
1237  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed || mkey.GetNVarCoeff())
1238  {
1239  NekDouble one = 1.0;
1240  DNekMatSharedPtr mat = GenMatrix(mkey);
1241 
1242  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,mat);
1243  }
1244  else
1245  {
1246  NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
1247  Array<TwoD, const NekDouble> df = m_metricinfo->GetDerivFactors(ptsKeys);
1248  int dir = 0;
1249  switch(mkey.GetMatrixType())
1250  {
1252  dir = 0;
1253  break;
1255  dir = 1;
1256  break;
1258  dir = 2;
1259  break;
1260  default:
1261  break;
1262  }
1263 
1265  mkey.GetShapeType(), *this);
1267  mkey.GetShapeType(), *this);
1268 
1269  DNekMat &deriv0 = *GetStdMatrix(deriv0key);
1270  DNekMat &deriv1 = *GetStdMatrix(deriv1key);
1271 
1272  int rows = deriv0.GetRows();
1273  int cols = deriv1.GetColumns();
1274 
1276  (*WeakDeriv) = df[2*dir][0]*deriv0 + df[2*dir+1][0]*deriv1;
1277 
1278  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(jac,WeakDeriv);
1279  }
1280  }
1281  break;
1283  {
1284  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed ||
1285  mkey.GetNVarCoeff())
1286  {
1287  NekDouble one = 1.0;
1288  DNekMatSharedPtr mat = GenMatrix(mkey);
1289 
1290  returnval = MemoryManager<DNekScalMat>::
1291  AllocateSharedPtr(one,mat);
1292  }
1293  else
1294  {
1295  int shapedim = 2;
1296 
1297  // dfdirxi = tan_{xi_x} * d \xi/dx
1298  // + tan_{xi_y} * d \xi/dy
1299  // + tan_{xi_z} * d \xi/dz
1300  // dfdireta = tan_{eta_x} * d \eta/dx
1301  // + tan_{xi_y} * d \xi/dy
1302  // + tan_{xi_z} * d \xi/dz
1303  NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
1305  m_metricinfo->GetDerivFactors(ptsKeys);
1306 
1307  Array<OneD, NekDouble> direction =
1309 
1310  // d / dx = df[0]*deriv0 + df[1]*deriv1
1311  // d / dy = df[2]*deriv0 + df[3]*deriv1
1312  // d / dz = df[4]*deriv0 + df[5]*deriv1
1313 
1314  // dfdir[dir] = e \cdot (d/dx, d/dy, d/dz)
1315  // = (e^0 * df[0] + e^1 * df[2]
1316  // + e^2 * df[4]) * deriv0
1317  // + (e^0 * df[1] + e^1 * df[3]
1318  // + e^2 * df[5]) * deriv1
1319  // dfdir[dir] = e^0 * df[2 * 0 + dir]
1320  // + e^1 * df[2 * 1 + dir]
1321  // + e^2 * df [ 2 * 2 + dir]
1322  Array<OneD, Array<OneD, NekDouble> > dfdir(shapedim);
1323  Expansion::ComputeGmatcdotMF(df, direction, dfdir);
1324 
1325  StdRegions::VarCoeffMap dfdirxi;
1326  StdRegions::VarCoeffMap dfdireta;
1327 
1328  dfdirxi[StdRegions::eVarCoeffWeakDeriv] = dfdir[0];
1329  dfdireta[StdRegions::eVarCoeffWeakDeriv] = dfdir[1];
1330 
1331  MatrixKey derivxikey( StdRegions::eWeakDeriv0,
1332  mkey.GetShapeType(), *this,
1334  dfdirxi);
1335  MatrixKey derivetakey( StdRegions::eWeakDeriv1,
1336  mkey.GetShapeType(), *this,
1338  dfdireta);
1339 
1340  DNekMat &derivxi = *GetStdMatrix(derivxikey);
1341  DNekMat &deriveta = *GetStdMatrix(derivetakey);
1342 
1343  int rows = derivxi.GetRows();
1344  int cols = deriveta.GetColumns();
1345 
1347  AllocateSharedPtr(rows,cols);
1348 
1349  (*WeakDirDeriv) = derivxi + deriveta;
1350 
1351  // Add (\nabla \cdot e^k ) Mass
1352  StdRegions::VarCoeffMap DiveMass;
1353  DiveMass[StdRegions::eVarCoeffMass] =
1355  StdRegions::StdMatrixKey stdmasskey(
1357  mkey.GetShapeType(), *this,
1359  DiveMass);
1360 
1361  DNekMatSharedPtr DiveMassmat = GetStdMatrix(stdmasskey);
1362 
1363  (*WeakDirDeriv) = (*WeakDirDeriv) + (*DiveMassmat);
1364 
1365  returnval = MemoryManager<DNekScalMat>::
1366  AllocateSharedPtr(jac,WeakDirDeriv);
1367  }
1368  break;
1369  }
1371  {
1372  if( (m_metricinfo->GetGtype() == SpatialDomains::eDeformed) ||
1374  {
1375  NekDouble one = 1.0;
1376  DNekMatSharedPtr mat = GenMatrix(mkey);
1377 
1378  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,mat);
1379  }
1380  else
1381  {
1383  mkey.GetShapeType(), *this);
1385  mkey.GetShapeType(), *this);
1387  mkey.GetShapeType(), *this);
1388 
1389  DNekMat &lap00 = *GetStdMatrix(lap00key);
1390  DNekMat &lap01 = *GetStdMatrix(lap01key);
1391  DNekMat &lap11 = *GetStdMatrix(lap11key);
1392 
1393  NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
1395  m_metricinfo->GetGmat(ptsKeys);
1396 
1397  int rows = lap00.GetRows();
1398  int cols = lap00.GetColumns();
1399 
1401 
1402  (*lap) = gmat[0][0] * lap00 +
1403  gmat[1][0] * (lap01 + Transpose(lap01)) +
1404  gmat[3][0] * lap11;
1405 
1406  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(jac,lap);
1407  }
1408  }
1409  break;
1411  {
1412  DNekMatSharedPtr mat = GenMatrix(mkey);
1413  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(1.0,mat);
1414  }
1415  break;
1417  {
1419 
1420  MatrixKey masskey(mkey, StdRegions::eMass);
1421  DNekScalMat &MassMat = *(this->m_matrixManager[masskey]);
1422 
1423  MatrixKey lapkey(mkey, StdRegions::eLaplacian);
1424  DNekScalMat &LapMat = *(this->m_matrixManager[lapkey]);
1425 
1426  int rows = LapMat.GetRows();
1427  int cols = LapMat.GetColumns();
1428 
1430 
1431  NekDouble one = 1.0;
1432  (*helm) = LapMat + factor*MassMat;
1433 
1434  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,helm);
1435  }
1436  break;
1438  {
1439  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1440  {
1441  NekDouble one = 1.0;
1442  DNekMatSharedPtr mat = GenMatrix(mkey);
1443  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,mat);
1444  }
1445  else
1446  {
1447  NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
1448  DNekMatSharedPtr mat = GetStdMatrix(mkey);
1449  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(jac,mat);
1450  }
1451  }
1452  break;
1456  {
1457  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1458  {
1459  NekDouble one = 1.0;
1460  DNekMatSharedPtr mat = GenMatrix(mkey);
1461  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,mat);
1462  }
1463  else
1464  {
1465  NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
1466 
1467  const Array<TwoD, const NekDouble>& df = m_metricinfo->GetDerivFactors(ptsKeys);
1468  int dir = 0;
1469 
1470  switch(mkey.GetMatrixType())
1471  {
1473  dir = 0;
1474  break;
1476  dir = 1;
1477  break;
1479  dir = 2;
1480  break;
1481  default:
1482  break;
1483  }
1484 
1486  mkey.GetShapeType(), *this);
1488  mkey.GetShapeType(), *this);
1489 
1490  DNekMat &stdiprod0 = *GetStdMatrix(iProdDeriv0Key);
1491  DNekMat &stdiprod1 = *GetStdMatrix(iProdDeriv0Key);
1492 
1493  int rows = stdiprod0.GetRows();
1494  int cols = stdiprod1.GetColumns();
1495 
1497  (*mat) = df[2*dir][0]*stdiprod0 + df[2*dir+1][0]*stdiprod1;
1498 
1499  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(jac,mat);
1500  }
1501  }
1502  break;
1503 
1505  {
1506  NekDouble one = 1.0;
1507 
1509 
1510  DNekMatSharedPtr mat = GenMatrix(hkey);
1511 
1512  mat->Invert();
1513  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,mat);
1514  }
1515  break;
1517  {
1518  NekDouble one = 1.0;
1519  MatrixKey helmkey(StdRegions::eHelmholtz, mkey.GetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
1520  DNekScalBlkMatSharedPtr helmStatCond = GetLocStaticCondMatrix(helmkey);
1521  DNekScalMatSharedPtr A =helmStatCond->GetBlock(0,0);
1523 
1525  }
1526  break;
1527  default:
1528  {
1529  NekDouble one = 1.0;
1530  DNekMatSharedPtr mat = GenMatrix(mkey);
1531 
1532  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,mat);
1533  }
1534  break;
1535  }
1536 
1537  return returnval;
1538  }
1539 
1540 
1542  {
1543  DNekScalBlkMatSharedPtr returnval;
1545 
1546  ASSERTL2(m_metricinfo->GetGtype() != SpatialDomains::eNoGeomType,"Geometric information is not set up");
1547 
1548  // set up block matrix system
1549  unsigned int nbdry = NumBndryCoeffs();
1550  unsigned int nint = (unsigned int)(m_ncoeffs - nbdry);
1551  unsigned int exp_size[] = {nbdry,nint};
1552  unsigned int nblks = 2;
1553  returnval = MemoryManager<DNekScalBlkMat>::AllocateSharedPtr(nblks,nblks,exp_size,exp_size); //Really need a constructor which takes Arrays
1554  NekDouble factor = 1.0;
1555 
1556  switch(mkey.GetMatrixType())
1557  {
1558  // this can only use stdregions statically condensed system for mass matrix
1559  case StdRegions::eMass:
1560  if((m_metricinfo->GetGtype() == SpatialDomains::eDeformed)||(mkey.GetNVarCoeff()))
1561  {
1562  factor = 1.0;
1563  goto UseLocRegionsMatrix;
1564  }
1565  else
1566  {
1567  factor = (m_metricinfo->GetJac(ptsKeys))[0];
1568  goto UseStdRegionsMatrix;
1569  }
1570  break;
1571  default: // use Deformed case for both regular and deformed geometries
1572  factor = 1.0;
1573  goto UseLocRegionsMatrix;
1574  break;
1575  UseStdRegionsMatrix:
1576  {
1577  NekDouble invfactor = 1.0/factor;
1578  NekDouble one = 1.0;
1580  DNekScalMatSharedPtr Atmp;
1581  DNekMatSharedPtr Asubmat;
1582 
1583  returnval->SetBlock(0,0,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(factor,Asubmat = mat->GetBlock(0,0)));
1584  returnval->SetBlock(0,1,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,Asubmat = mat->GetBlock(0,1)));
1585  returnval->SetBlock(1,0,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(factor,Asubmat = mat->GetBlock(1,0)));
1586  returnval->SetBlock(1,1,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(invfactor,Asubmat = mat->GetBlock(1,1)));
1587  }
1588  break;
1589 
1590  UseLocRegionsMatrix:
1591  {
1592  int i,j;
1593  NekDouble invfactor = 1.0/factor;
1594  NekDouble one = 1.0;
1595 
1596  DNekScalMat &mat = *GetLocMatrix(mkey);
1597 
1602 
1603  Array<OneD,unsigned int> bmap(nbdry);
1604  Array<OneD,unsigned int> imap(nint);
1605  GetBoundaryMap(bmap);
1606  GetInteriorMap(imap);
1607 
1608  for(i = 0; i < nbdry; ++i)
1609  {
1610  for(j = 0; j < nbdry; ++j)
1611  {
1612  (*A)(i,j) = mat(bmap[i],bmap[j]);
1613  }
1614 
1615  for(j = 0; j < nint; ++j)
1616  {
1617  (*B)(i,j) = mat(bmap[i],imap[j]);
1618  }
1619  }
1620 
1621  for(i = 0; i < nint; ++i)
1622  {
1623  for(j = 0; j < nbdry; ++j)
1624  {
1625  (*C)(i,j) = mat(imap[i],bmap[j]);
1626  }
1627 
1628  for(j = 0; j < nint; ++j)
1629  {
1630  (*D)(i,j) = mat(imap[i],imap[j]);
1631  }
1632  }
1633 
1634  // Calculate static condensed system
1635  if(nint)
1636  {
1637  D->Invert();
1638  (*B) = (*B)*(*D);
1639  (*A) = (*A) - (*B)*(*C);
1640  }
1641 
1642  DNekScalMatSharedPtr Atmp;
1643 
1644  returnval->SetBlock(0,0,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(factor,A));
1645  returnval->SetBlock(0,1,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,B));
1646  returnval->SetBlock(1,0,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(factor,C));
1647  returnval->SetBlock(1,1,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(invfactor,D));
1648 
1649  }
1650  }
1651 
1652  return returnval;
1653  }
1654 
1655 
1657  {
1658  return m_matrixManager[mkey];
1659  }
1660 
1661 
1663  {
1664  return m_staticCondMatrixManager[mkey];
1665  }
1666 
1668  {
1669  m_staticCondMatrixManager.DeleteObject(mkey);
1670  }
1671 
1672 
1673 
1675  Array<OneD,NekDouble> &outarray,
1676  const StdRegions::StdMatrixKey &mkey)
1677  {
1678  StdExpansion::MassMatrixOp_MatFree(inarray,outarray,mkey);
1679  }
1680 
1681 
1683  Array<OneD,NekDouble> &outarray,
1684  const StdRegions::StdMatrixKey &mkey)
1685  {
1686  TriExp::LaplacianMatrixOp_MatFree(inarray,outarray,mkey);
1687  }
1688 
1689 
1690  void TriExp::v_LaplacianMatrixOp(const int k1, const int k2,
1691  const Array<OneD, const NekDouble> &inarray,
1692  Array<OneD,NekDouble> &outarray,
1693  const StdRegions::StdMatrixKey &mkey)
1694  {
1695  StdExpansion::LaplacianMatrixOp_MatFree(k1,k2,inarray,outarray,mkey);
1696  }
1697 
1698 
1700  const Array<OneD, const NekDouble> &inarray,
1701  Array<OneD,NekDouble> &outarray,
1702  const StdRegions::StdMatrixKey &mkey)
1703  {
1704  StdExpansion::WeakDerivMatrixOp_MatFree(i,inarray,outarray,mkey);
1705  }
1706 
1707 
1709  Array<OneD,NekDouble> &outarray,
1710  const StdRegions::StdMatrixKey &mkey)
1711  {
1712  StdExpansion::WeakDirectionalDerivMatrixOp_MatFree(inarray,outarray,mkey);
1713  }
1714 
1715 
1717  Array<OneD,NekDouble> &outarray,
1718  const StdRegions::StdMatrixKey &mkey)
1719  {
1720  StdExpansion::MassLevelCurvatureMatrixOp_MatFree(inarray,outarray,mkey);
1721  }
1722 
1723 
1725  Array<OneD,NekDouble> &outarray,
1726  const StdRegions::StdMatrixKey &mkey)
1727  {
1728  TriExp::HelmholtzMatrixOp_MatFree(inarray,outarray,mkey);
1729  }
1730 
1731 
1733  Array<OneD,NekDouble> &outarray,
1734  const StdRegions::StdMatrixKey &mkey)
1735  {
1736  DNekScalMatSharedPtr mat = GetLocMatrix(mkey);
1737 
1738  if(inarray.get() == outarray.get())
1739  {
1741  Vmath::Vcopy(m_ncoeffs,inarray.get(),1,tmp.get(),1);
1742 
1743  Blas::Dgemv('N',m_ncoeffs,m_ncoeffs,mat->Scale(),(mat->GetOwnedMatrix())->GetPtr().get(),
1744  m_ncoeffs, tmp.get(), 1, 0.0, outarray.get(), 1);
1745  }
1746  else
1747  {
1748  Blas::Dgemv('N',m_ncoeffs,m_ncoeffs,mat->Scale(),(mat->GetOwnedMatrix())->GetPtr().get(),
1749  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
1750  }
1751  }
1752 
1753 
1755  const Array<OneD, const NekDouble> &inarray,
1756  Array<OneD, NekDouble> &outarray,
1758  {
1759  if (m_metrics.count(eMetricLaplacian00) == 0)
1760  {
1762  }
1763 
1764  int nquad0 = m_base[0]->GetNumPoints();
1765  int nquad1 = m_base[1]->GetNumPoints();
1766  int nqtot = nquad0*nquad1;
1767  int nmodes0 = m_base[0]->GetNumModes();
1768  int nmodes1 = m_base[1]->GetNumModes();
1769  int wspsize = max(max(max(nqtot,m_ncoeffs),nquad1*nmodes0),nquad0*nmodes1);
1770 
1771  ASSERTL1(wsp.size() >= 3*wspsize,
1772  "Workspace is of insufficient size.");
1773 
1774  const Array<OneD, const NekDouble>& base0 = m_base[0]->GetBdata();
1775  const Array<OneD, const NekDouble>& base1 = m_base[1]->GetBdata();
1776  const Array<OneD, const NekDouble>& dbase0 = m_base[0]->GetDbdata();
1777  const Array<OneD, const NekDouble>& dbase1 = m_base[1]->GetDbdata();
1781 
1782  // Allocate temporary storage
1783  Array<OneD,NekDouble> wsp0(wsp);
1784  Array<OneD,NekDouble> wsp1(wsp+wspsize);
1785  Array<OneD,NekDouble> wsp2(wsp+2*wspsize);
1786 
1787  StdExpansion2D::PhysTensorDeriv(inarray,wsp1,wsp2);
1788 
1789  // wsp0 = k = g0 * wsp1 + g1 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
1790  // wsp2 = l = g1 * wsp1 + g2 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
1791  // where g0, g1 and g2 are the metric terms set up in the GeomFactors class
1792  // especially for this purpose
1793  Vmath::Vvtvvtp(nqtot,&metric00[0],1,&wsp1[0],1,&metric01[0],1,&wsp2[0],1,&wsp0[0],1);
1794  Vmath::Vvtvvtp(nqtot,&metric01[0],1,&wsp1[0],1,&metric11[0],1,&wsp2[0],1,&wsp2[0],1);
1795 
1796  // outarray = m = (D_xi1 * B)^T * k
1797  // wsp1 = n = (D_xi2 * B)^T * l
1798  IProductWRTBase_SumFacKernel(dbase0,base1,wsp0,outarray,wsp1);
1799  IProductWRTBase_SumFacKernel(base0,dbase1,wsp2,wsp1, wsp0);
1800 
1801  // outarray = outarray + wsp1
1802  // = L * u_hat
1803  Vmath::Vadd(m_ncoeffs,wsp1.get(),1,outarray.get(),1,outarray.get(),1);
1804  }
1805 
1806 
1808  {
1809  if (m_metrics.count(eMetricQuadrature) == 0)
1810  {
1812  }
1813 
1814  unsigned int i, j;
1815  const SpatialDomains::GeomType type = m_metricinfo->GetGtype();
1816  const unsigned int nqtot = GetTotPoints();
1817  const unsigned int dim = 2;
1821  };
1822 
1823  Array<OneD, NekDouble> dEta_dXi[2] = {Array<OneD, NekDouble>(nqtot,1.0),
1824  Array<OneD, NekDouble>(nqtot,1.0)};
1825 
1826  for (i = 0; i < dim; ++i)
1827  {
1828  for (j = i; j < dim; ++j)
1829  {
1830  m_metrics[m[i][j]] = Array<OneD, NekDouble>(nqtot);
1831  }
1832  }
1833 
1834  const Array<OneD, const NekDouble>& z0 = m_base[0]->GetZ();
1835  const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
1836  const unsigned int nquad0 = m_base[0]->GetNumPoints();
1837  const unsigned int nquad1 = m_base[1]->GetNumPoints();
1838  const Array<TwoD, const NekDouble>& df =
1839  m_metricinfo->GetDerivFactors(GetPointsKeys());
1840 
1841  for(i = 0; i < nquad1; i++)
1842  {
1843  Blas::Dscal(nquad0,2.0/(1-z1[i]),&dEta_dXi[0][0]+i*nquad0,1);
1844  Blas::Dscal(nquad0,2.0/(1-z1[i]),&dEta_dXi[1][0]+i*nquad0,1);
1845  }
1846  for(i = 0; i < nquad0; i++)
1847  {
1848  Blas::Dscal(nquad1,0.5*(1+z0[i]),&dEta_dXi[1][0]+i,nquad0);
1849  }
1850 
1851  Array<OneD, NekDouble> tmp(nqtot);
1852  if((type == SpatialDomains::eRegular ||
1854  {
1855  Vmath::Smul (nqtot,df[0][0],&dEta_dXi[0][0],1,&tmp[0],1);
1856  Vmath::Svtvp(nqtot,df[1][0],&dEta_dXi[1][0],1,&tmp[0],1,&tmp[0],1);
1857 
1858  Vmath::Vmul (nqtot,&tmp[0],1, &tmp[0],1,&m_metrics[eMetricLaplacian00][0],1);
1859  Vmath::Smul (nqtot,df[1][0],&tmp[0],1,&m_metrics[eMetricLaplacian01][0],1);
1860 
1861 
1862  Vmath::Smul (nqtot,df[2][0],&dEta_dXi[0][0],1,&tmp[0],1);
1863  Vmath::Svtvp(nqtot,df[3][0],&dEta_dXi[1][0],1,&tmp[0],1,&tmp[0],1);
1864 
1865  Vmath::Vvtvp(nqtot,&tmp[0],1, &tmp[0],1,&m_metrics[eMetricLaplacian00][0],1,&m_metrics[eMetricLaplacian00][0],1);
1866  Vmath::Svtvp(nqtot,df[3][0],&tmp[0],1,&m_metrics[eMetricLaplacian01][0],1,&m_metrics[eMetricLaplacian01][0],1);
1867 
1868  if(GetCoordim() == 3)
1869  {
1870  Vmath::Smul (nqtot,df[4][0],&dEta_dXi[0][0],1,&tmp[0],1);
1871  Vmath::Svtvp(nqtot,df[5][0],&dEta_dXi[1][0],1,&tmp[0],1,&tmp[0],1);
1872 
1873  Vmath::Vvtvp(nqtot,&tmp[0],1, &tmp[0],1,&m_metrics[eMetricLaplacian00][0],1,&m_metrics[eMetricLaplacian00][0],1);
1874  Vmath::Svtvp(nqtot,df[5][0],&tmp[0],1,&m_metrics[eMetricLaplacian01][0],1,&m_metrics[eMetricLaplacian01][0],1);
1875  }
1876 
1877  NekDouble g2 = df[1][0]*df[1][0] + df[3][0]*df[3][0];
1878  if(GetCoordim() == 3)
1879  {
1880  g2 += df[5][0]*df[5][0];
1881  }
1882  Vmath::Fill(nqtot,g2,&m_metrics[eMetricLaplacian11][0],1);
1883  }
1884  else
1885  {
1886 
1887  Vmath::Vmul (nqtot,&df[0][0],1,&dEta_dXi[0][0],1,&tmp[0],1);
1888  Vmath::Vvtvp(nqtot,&df[1][0],1,&dEta_dXi[1][0],1,&tmp[0],1,&tmp[0],1);
1889 
1890  Vmath::Vmul (nqtot,&tmp[0], 1,&tmp[0], 1,&m_metrics[eMetricLaplacian00][0],1);
1891  Vmath::Vmul (nqtot,&df[1][0],1,&tmp[0], 1,&m_metrics[eMetricLaplacian01][0],1);
1892  Vmath::Vmul (nqtot,&df[1][0],1,&df[1][0],1,&m_metrics[eMetricLaplacian11][0],1);
1893 
1894 
1895  Vmath::Vmul (nqtot,&df[2][0],1,&dEta_dXi[0][0],1,&tmp[0],1);
1896  Vmath::Vvtvp(nqtot,&df[3][0],1,&dEta_dXi[1][0],1,&tmp[0],1,&tmp[0],1);
1897 
1898  Vmath::Vvtvp(nqtot,&tmp[0], 1,&tmp[0], 1,&m_metrics[eMetricLaplacian00][0],1,&m_metrics[eMetricLaplacian00][0],1);
1899  Vmath::Vvtvp(nqtot,&df[3][0],1,&tmp[0], 1,&m_metrics[eMetricLaplacian01][0],1,&m_metrics[eMetricLaplacian01][0],1);
1900  Vmath::Vvtvp(nqtot,&df[3][0],1,&df[3][0],1,&m_metrics[eMetricLaplacian11][0],1,&m_metrics[eMetricLaplacian11][0],1);
1901 
1902  if(GetCoordim() == 3)
1903  {
1904  Vmath::Vmul (nqtot,&df[4][0],1,&dEta_dXi[0][0],1,&tmp[0],1);
1905  Vmath::Vvtvp(nqtot,&df[5][0],1,&dEta_dXi[1][0],1,&tmp[0],1,&tmp[0],1);
1906 
1907  Vmath::Vvtvp(nqtot,&tmp[0], 1,&tmp[0], 1,&m_metrics[eMetricLaplacian00][0],1,&m_metrics[eMetricLaplacian00][0],1);
1908  Vmath::Vvtvp(nqtot,&df[5][0],1,&tmp[0], 1,&m_metrics[eMetricLaplacian01][0],1,&m_metrics[eMetricLaplacian01][0],1);
1909  Vmath::Vvtvp(nqtot,&df[5][0],1,&df[5][0],1,&m_metrics[eMetricLaplacian11][0],1,&m_metrics[eMetricLaplacian11][0],1);
1910  }
1911  }
1912 
1913  for (unsigned int i = 0; i < dim; ++i)
1914  {
1915  for (unsigned int j = i; j < dim; ++j)
1916  {
1918  m_metrics[m[i][j]]);
1919 
1920  }
1921  }
1922  }
1923 
1924  /**
1925  * Function is used to compute exactly the advective numerical flux on
1926  * theinterface of two elements with different expansions, hence an
1927  * appropriate number of Gauss points has to be used. The number of
1928  * Gauss points has to be equal to the number used by the highest
1929  * polynomial degree of the two adjacent elements. Furthermore, this
1930  * function is used to compute the sensor value in each element.
1931  *
1932  * @param numMin Is the reduced polynomial order
1933  * @param inarray Input array of coefficients
1934  * @param dumpVar Output array of reduced coefficients.
1935  */
1937  int numMin,
1938  const Array<OneD, const NekDouble> &inarray,
1939  Array<OneD, NekDouble> &outarray)
1940  {
1941  int n_coeffs = inarray.size();
1942  int nquad0 = m_base[0]->GetNumPoints();
1943  int nquad1 = m_base[1]->GetNumPoints();
1944  int nqtot = nquad0*nquad1;
1945  int nmodes0 = m_base[0]->GetNumModes();
1946  int nmodes1 = m_base[1]->GetNumModes();
1947  int numMin2 = nmodes0, i;
1948 
1949  Array<OneD, NekDouble> coeff(n_coeffs,0.0);
1950  Array<OneD, NekDouble> phys_tmp(nqtot,0.0);
1951  Array<OneD, NekDouble> tmp, tmp2;
1952 
1953  const LibUtilities::PointsKey Pkey0 = m_base[0]->GetPointsKey();
1954  const LibUtilities::PointsKey Pkey1 = m_base[1]->GetPointsKey();
1955 
1957  m_base[0]->GetBasisType(), nmodes0, Pkey0);
1959  m_base[1]->GetBasisType(), nmodes1, Pkey1);
1960  LibUtilities::BasisKey bortho0(
1961  LibUtilities::eOrtho_A, nmodes0, Pkey0);
1962  LibUtilities::BasisKey bortho1(
1963  LibUtilities::eOrtho_B, nmodes1, Pkey1);
1964 
1965  // Check if it is also possible to use the same InterCoeff routine
1966  // which is also used for Quadrilateral and Hexagonal shaped
1967  // elements
1968 
1969  // For now, set up the used basis on the standard element to
1970  // calculate the phys values, set up the orthogonal basis to do a
1971  // forward transform, to obtain the coefficients in orthogonal
1972  // coefficient space
1973  StdRegions::StdTriExpSharedPtr m_OrthoTriExp;
1975 
1977  ::AllocateSharedPtr(b0, b1);
1978  m_OrthoTriExp = MemoryManager<StdRegions::StdTriExp>
1979  ::AllocateSharedPtr(bortho0, bortho1);
1980 
1981  m_TriExp ->BwdTrans(inarray,phys_tmp);
1982  m_OrthoTriExp->FwdTrans(phys_tmp, coeff);
1983 
1984  for (i = 0; i < n_coeffs; i++)
1985  {
1986  if (i == numMin)
1987  {
1988  coeff[i] = 0.0;
1989  numMin += numMin2 - 1;
1990  numMin2 -= 1.0;
1991  }
1992  }
1993 
1994  m_OrthoTriExp->BwdTrans(coeff,phys_tmp);
1995  m_TriExp ->FwdTrans(phys_tmp, outarray);
1996  }
1997 
1999  Array<OneD, NekDouble> &array,
2000  const StdRegions::StdMatrixKey &mkey)
2001  {
2002  int nq = GetTotPoints();
2003 
2004  // Calculate sqrt of the Jacobian
2006  m_metricinfo->GetJac(GetPointsKeys());
2007  Array<OneD, NekDouble> sqrt_jac(nq);
2008  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
2009  {
2010  Vmath::Vsqrt(nq,jac,1,sqrt_jac,1);
2011  }
2012  else
2013  {
2014  Vmath::Fill(nq,sqrt(jac[0]),sqrt_jac,1);
2015  }
2016 
2017  // Multiply array by sqrt(Jac)
2018  Vmath::Vmul(nq,sqrt_jac,1,array,1,array,1);
2019 
2020  // Apply std region filter
2021  StdTriExp::v_SVVLaplacianFilter( array, mkey);
2022 
2023  // Divide by sqrt(Jac)
2024  Vmath::Vdiv(nq,array,1,sqrt_jac,1,array,1);
2025  }
2026 
2027  }
2028 }
2029 
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:274
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:15
Describes the specification for a Basis.
Definition: Basis.h:50
Defines a specification for a set of points.
Definition: Points.h:60
unsigned int GetNumPoints() const
Definition: Points.h:107
virtual DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey)
SpatialDomains::Geometry2DSharedPtr GetGeom2D() const
Definition: Expansion2D.h:204
std::map< int, NormalVector > m_edgeNormals
Definition: Expansion2D.h:121
std::map< int, Array< OneD, NekDouble > > m_elmtBndNormDirElmtLen
the element length in each element boundary(Vertex, edge or face) normal direction calculated based o...
Definition: Expansion.h:284
SpatialDomains::GeometrySharedPtr GetGeom() const
Definition: Expansion.cpp:172
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:272
ExpansionSharedPtr GetLeftAdjacentElementExp() const
Definition: Expansion.h:438
DNekMatSharedPtr BuildVertexMatrix(const DNekScalMatSharedPtr &r_bnd)
Definition: Expansion.cpp:103
void ComputeGmatcdotMF(const Array< TwoD, const NekDouble > &df, const Array< OneD, const NekDouble > &direction, Array< OneD, Array< OneD, NekDouble > > &dfdir)
Definition: Expansion.cpp:399
int GetLeftAdjacentElementTrace() const
Definition: Expansion.h:454
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:273
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3)
Definition: Expansion.cpp:318
StdRegions::Orientation GetTraceOrient(int trace)
Definition: Expansion.h:164
DNekScalMatSharedPtr GetLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:90
virtual void v_ComputeTraceNormal(const int edge)
Definition: TriExp.cpp:909
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: TriExp.h:291
virtual void v_PhysDirectionalDeriv(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &out)
Physical derivative along a direction vector.
Definition: TriExp.cpp:196
virtual void v_NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
Definition: TriExp.cpp:684
virtual void v_GetTracePhysMap(const int edge, Array< OneD, int > &outarray)
Definition: TriExp.cpp:870
virtual void v_GetCoord(const Array< OneD, const NekDouble > &Lcoords, Array< OneD, NekDouble > &coords)
Definition: TriExp.cpp:739
virtual void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
Definition: TriExp.cpp:1682
virtual void v_ExtractDataToCoeffs(const NekDouble *data, const std::vector< unsigned int > &nummodes, const int mode_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType)
Definition: TriExp.cpp:1091
virtual void v_WeakDerivMatrixOp(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
Definition: TriExp.cpp:1699
virtual void v_MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
Definition: TriExp.cpp:1674
virtual DNekScalMatSharedPtr v_GetLocMatrix(const MatrixKey &mkey)
Definition: TriExp.cpp:1656
virtual void v_ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: TriExp.cpp:1936
virtual const LibUtilities::BasisSharedPtr & v_GetBasis(int dir) const
Definition: TriExp.cpp:1139
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: TriExp.cpp:414
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix(const MatrixKey &mkey)
Definition: TriExp.cpp:1662
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: TriExp.cpp:108
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: TriExp.h:292
virtual int v_GetCoordim()
Definition: TriExp.cpp:1085
virtual int v_GetNumPoints(const int dir) const
Definition: TriExp.cpp:1146
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
Definition: TriExp.cpp:422
virtual void v_GetEdgeInterpVals(const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: TriExp.cpp:850
virtual StdRegions::Orientation v_GetTraceOrient(int edge)
Definition: TriExp.cpp:1134
virtual StdRegions::StdExpansionSharedPtr v_GetStdExp(void) const
Definition: TriExp.cpp:720
virtual void v_GetTraceQFactors(const int edge, Array< OneD, NekDouble > &outarray)
Definition: TriExp.cpp:859
void v_DropLocStaticCondMatrix(const MatrixKey &mkey)
Definition: TriExp.cpp:1667
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray)
Integrates the specified function over the domain.
Definition: TriExp.cpp:85
virtual void v_MassLevelCurvatureMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
Definition: TriExp.cpp:1716
virtual void v_ComputeLaplacianMetric()
Definition: TriExp.cpp:1807
virtual DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey)
Definition: TriExp.cpp:1152
virtual void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdRegions::StdMatrixKey &mkey)
Definition: TriExp.cpp:1998
virtual void v_LaplacianMatrixOp_MatFree_Kernel(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
Definition: TriExp.cpp:1754
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdRegions::StdMatrixKey &mkey)
Definition: TriExp.cpp:1175
virtual void v_IProductWRTBase_MatOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: TriExp.cpp:448
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: TriExp.cpp:267
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: TriExp.cpp:407
virtual NekDouble v_PhysEvaluate(const Array< OneD, const NekDouble > &coord, const Array< OneD, const NekDouble > &physvals)
This function evaluates the expansion at a single (arbitrary) point of the domain.
Definition: TriExp.cpp:778
virtual NekDouble v_StdPhysEvaluate(const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
Definition: TriExp.cpp:770
virtual DNekScalMatSharedPtr CreateMatrix(const MatrixKey &mkey)
Definition: TriExp.cpp:1186
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3)
Definition: TriExp.cpp:756
virtual void v_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
Definition: TriExp.cpp:1724
virtual void v_GetTracePhysVals(const int edge, const StdRegions::StdExpansionSharedPtr &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient)
Definition: TriExp.cpp:789
virtual void v_IProductWRTDirectionalDerivBase(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: TriExp.cpp:600
virtual void v_FwdTrans_BndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: TriExp.cpp:285
virtual StdRegions::StdExpansionSharedPtr v_GetLinStdExp(void) const
Definition: TriExp.cpp:728
virtual void v_IProductWRTDirectionalDerivBase_SumFac(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Directinoal Derivative in the modal space in the dir direction of varcoeffs.
Definition: TriExp.cpp:614
virtual DNekScalBlkMatSharedPtr CreateStaticCondMatrix(const MatrixKey &mkey)
Definition: TriExp.cpp:1541
virtual void v_GeneralMatrixOp_MatOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
Definition: TriExp.cpp:1732
virtual void v_AlignVectorToCollapsedDir(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: TriExp.cpp:494
virtual void v_WeakDirectionalDerivMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
Definition: TriExp.cpp:1708
virtual void v_IProductWRTDerivBase_MatOp(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: TriExp.cpp:560
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: TriExp.cpp:461
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
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)
void GetBoundaryMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:687
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:134
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:158
void HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:617
const LibUtilities::PointsKeyVector GetPointsKeys() const
DNekBlkMatSharedPtr GetStdStaticCondMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:622
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:762
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix(const LocalRegions::MatrixKey &mkey)
Definition: StdExpansion.h:660
void IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:628
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:537
void LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void GetTraceToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
Definition: StdExpansion.h:703
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:376
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:692
DNekMatSharedPtr GenMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:850
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:208
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:221
void IProductWRTDirectionalDerivBase_SumFac(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:733
void PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Definition: StdExpansion.h:855
Array< OneD, LibUtilities::BasisSharedPtr > m_base
const ConstFactorMap & GetConstFactors() const
Definition: StdMatrixKey.h:135
LibUtilities::ShapeType GetShapeType() const
Definition: StdMatrixKey.h:86
const VarCoeffMap & GetVarCoeffs() const
Definition: StdMatrixKey.h:161
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:81
const Array< OneD, const NekDouble > & GetVarCoeff(const StdRegions::VarCoeffType &coeff) const
Definition: StdMatrixKey.h:145
NekDouble GetConstFactor(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:121
bool ConstFactorExists(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:130
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:265
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition: Blas.hpp:182
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:167
std::shared_ptr< Basis > BasisSharedPtr
void Interp1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
this function interpolates a 1D function evaluated at the quadrature points of the basis fbasis0 to ...
Definition: Interp.cpp:53
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:246
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:49
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:45
@ eOrtho_B
Principle Orthogonal Functions .
Definition: BasisType.h:46
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48
std::shared_ptr< SegExp > SegExpSharedPtr
Definition: SegExp.h:272
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition: GeomFactors.h:62
GeomType
Indicates the type of element geometry.
@ eRegular
Geometry is straight-sided with constant geometric factors.
@ eNoGeomType
No type defined.
@ eMovingRegular
Currently unused.
@ eDeformed
Geometry is curved or has non-constant factors.
std::shared_ptr< TriGeom > TriGeomSharedPtr
Definition: TriGeom.h:58
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:272
std::shared_ptr< StdTriExp > StdTriExpSharedPtr
Definition: StdTriExp.h:272
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:315
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:71
NekMatrix< InnerMatrixType, BlockMatrixTag > Transpose(NekMatrix< InnerMatrixType, BlockMatrixTag > &rhs)
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:73
static Array< OneD, NekDouble > NullNekDouble1DArray
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
double NekDouble
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
Definition: Vmath.cpp:475
void Svtsvtp(int n, const T alpha, const T *x, int incx, const T beta, const T *y, int incy, T *z, int incz)
vvtvvtp (scalar times vector plus scalar times vector):
Definition: Vmath.cpp:691
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:192
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
Definition: Vmath.cpp:565
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:513
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:322
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:225
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:291
void Vdiv(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:257
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:436
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45
void Reverse(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1226
void Vvtvvtp(int n, const T *v, int incv, const T *w, int incw, const T *x, int incx, const T *y, int incy, T *z, int incz)
vvtvvtp (vector times vector plus vector times vector):
Definition: Vmath.cpp:625
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199
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:372
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:267