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