Nektar++
StdPyrExp.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: StdPyrExp.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: pyramidic routines built upon StdExpansion3D
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <boost/core/ignore_unused.hpp>
36 
38 #include <StdRegions/StdPyrExp.h>
39 #include <iomanip>
40 
41 using namespace std;
42 
43 namespace Nektar
44 {
45 namespace StdRegions
46 {
47 StdPyrExp::StdPyrExp(const LibUtilities::BasisKey &Ba,
48  const LibUtilities::BasisKey &Bb,
49  const LibUtilities::BasisKey &Bc)
50  : StdExpansion(LibUtilities::StdPyrData::getNumberOfCoefficients(
51  Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
52  3, Ba, Bb, Bc),
53  StdExpansion3D(LibUtilities::StdPyrData::getNumberOfCoefficients(
54  Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
55  Ba, Bb, Bc)
56 {
57 
58  ASSERTL0(Ba.GetNumModes() <= Bc.GetNumModes(),
59  "order in 'a' direction is higher "
60  "than order in 'c' direction");
61  ASSERTL0(Bb.GetNumModes() <= Bc.GetNumModes(),
62  "order in 'b' direction is higher "
63  "than order in 'c' direction");
66  "Expected basis type in 'c' direction to be ModifiedPyr_C or "
67  "OrthoPyr_C");
68 }
69 
71 {
72 }
73 
74 //---------------------------------------
75 // Differentiation/integration Methods
76 //---------------------------------------
77 
78 /**
79  * \brief Calculate the derivative of the physical points
80  *
81  * The derivative is evaluated at the nodal physical points.
82  * Derivatives with respect to the local Cartesian coordinates.
83  *
84  * \f$\begin{Bmatrix} \frac {\partial} {\partial \xi_1} \\ \frac
85  * {\partial} {\partial \xi_2} \\ \frac {\partial} {\partial \xi_3}
86  * \end{Bmatrix} = \begin{Bmatrix} \frac 2 {(1-\eta_3)} \frac \partial
87  * {\partial \bar \eta_1} \\ \frac {\partial} {\partial \xi_2} \ \
88  * \frac {(1 + \bar \eta_1)} {(1 - \eta_3)} \frac \partial {\partial
89  * \bar \eta_1} + \frac {\partial} {\partial \eta_3} \end{Bmatrix}\f$
90  */
92  Array<OneD, NekDouble> &out_dxi1,
93  Array<OneD, NekDouble> &out_dxi2,
94  Array<OneD, NekDouble> &out_dxi3)
95 {
96  // PhysDerivative implementation based on Spen's book page 152.
97  int Qx = m_base[0]->GetNumPoints();
98  int Qy = m_base[1]->GetNumPoints();
99  int Qz = m_base[2]->GetNumPoints();
100 
101  Array<OneD, NekDouble> dEta_bar1(Qx * Qy * Qz, 0.0);
102  Array<OneD, NekDouble> dXi2(Qx * Qy * Qz, 0.0);
103  Array<OneD, NekDouble> dEta3(Qx * Qy * Qz, 0.0);
104  PhysTensorDeriv(u_physical, dEta_bar1, dXi2, dEta3);
105 
106  Array<OneD, const NekDouble> eta_x, eta_y, eta_z;
107  eta_x = m_base[0]->GetZ();
108  eta_y = m_base[1]->GetZ();
109  eta_z = m_base[2]->GetZ();
110 
111  int i, j, k, n;
112 
113  if (out_dxi1.size() > 0)
114  {
115  for (k = 0, n = 0; k < Qz; ++k)
116  {
117  NekDouble fac = 2.0 / (1.0 - eta_z[k]);
118  for (j = 0; j < Qy; ++j)
119  {
120  for (i = 0; i < Qx; ++i, ++n)
121  {
122  out_dxi1[n] = fac * dEta_bar1[n];
123  }
124  }
125  }
126  }
127 
128  if (out_dxi2.size() > 0)
129  {
130  for (k = 0, n = 0; k < Qz; ++k)
131  {
132  NekDouble fac = 2.0 / (1.0 - eta_z[k]);
133  for (j = 0; j < Qy; ++j)
134  {
135  for (i = 0; i < Qx; ++i, ++n)
136  {
137  out_dxi2[n] = fac * dXi2[n];
138  }
139  }
140  }
141  }
142 
143  if (out_dxi3.size() > 0)
144  {
145  for (k = 0, n = 0; k < Qz; ++k)
146  {
147  NekDouble fac = 1.0 / (1.0 - eta_z[k]);
148  for (j = 0; j < Qy; ++j)
149  {
150  NekDouble fac1 = (1.0 + eta_y[j]);
151  for (i = 0; i < Qx; ++i, ++n)
152  {
153  out_dxi3[n] = (1.0 + eta_x[i]) * fac * dEta_bar1[n] +
154  fac1 * fac * dXi2[n] + dEta3[n];
155  }
156  }
157  }
158  }
159 }
160 
161 void StdPyrExp::v_PhysDeriv(const int dir,
162  const Array<OneD, const NekDouble> &inarray,
163  Array<OneD, NekDouble> &outarray)
164 {
165  switch (dir)
166  {
167  case 0:
168  {
169  v_PhysDeriv(inarray, outarray, NullNekDouble1DArray,
171  break;
172  }
173 
174  case 1:
175  {
176  v_PhysDeriv(inarray, NullNekDouble1DArray, outarray,
178  break;
179  }
180 
181  case 2:
182  {
184  outarray);
185  break;
186  }
187 
188  default:
189  {
190  ASSERTL1(false, "input dir is out of range");
191  }
192  break;
193  }
194 }
195 
197  Array<OneD, NekDouble> &out_d0,
198  Array<OneD, NekDouble> &out_d1,
199  Array<OneD, NekDouble> &out_d2)
200 {
201  StdPyrExp::v_PhysDeriv(inarray, out_d0, out_d1, out_d2);
202 }
203 
204 void StdPyrExp::v_StdPhysDeriv(const int dir,
205  const Array<OneD, const NekDouble> &inarray,
206  Array<OneD, NekDouble> &outarray)
207 {
208  StdPyrExp::v_PhysDeriv(dir, inarray, outarray);
209 }
210 
211 //---------------------------------------
212 // Transforms
213 //---------------------------------------
214 
215 /**
216  * \brief Backward transformation is evaluated at the quadrature
217  * points.
218  *
219  * \f$ u^{\delta} (\xi_{1i}, \xi_{2j}, \xi_{3k}) = \sum_{m(pqr)} \hat
220  * u_{pqr} \phi_{pqr} (\xi_{1i}, \xi_{2j}, \xi_{3k})\f$
221  *
222  * Backward transformation is three dimensional tensorial expansion
223  *
224  * \f$ u (\xi_{1i}, \xi_{2j}, \xi_{3k}) = \sum_{p=0}^{Q_x} \psi_p^a
225  * (\xi_{1i}) \lbrace { \sum_{q=0}^{Q_y} \psi_{q}^a (\xi_{2j})
226  * \lbrace { \sum_{r=0}^{Q_z} \hat u_{pqr} \psi_{pqr}^c (\xi_{3k})
227  * \rbrace} \rbrace}. \f$
228  *
229  * And sumfactorizing step of the form is as:\ \ \f$ f_{pqr}
230  * (\xi_{3k}) = \sum_{r=0}^{Q_z} \hat u_{pqr} \psi_{pqr}^c
231  * (\xi_{3k}),\\ g_{p} (\xi_{2j}, \xi_{3k}) = \sum_{r=0}^{Q_y}
232  * \psi_{p}^a (\xi_{2j}) f_{pqr} (\xi_{3k}),\\ u(\xi_{1i}, \xi_{2j},
233  * \xi_{3k}) = \sum_{p=0}^{Q_x} \psi_{p}^a (\xi_{1i}) g_{p}
234  * (\xi_{2j}, \xi_{3k}). \f$
235  **/
237  Array<OneD, NekDouble> &outarray)
238 {
239  if (m_base[0]->Collocation() && m_base[1]->Collocation() &&
240  m_base[2]->Collocation())
241  {
243  m_base[2]->GetNumPoints(),
244  inarray, 1, outarray, 1);
245  }
246  else
247  {
248  StdPyrExp::v_BwdTrans_SumFac(inarray, outarray);
249  }
250 }
251 
252 /**
253  * Sum-factorisation implementation of the BwdTrans operation.
254  */
256  Array<OneD, NekDouble> &outarray)
257 {
258  int nquad0 = m_base[0]->GetNumPoints();
259  int nquad1 = m_base[1]->GetNumPoints();
260  int nquad2 = m_base[2]->GetNumPoints();
261  int order0 = m_base[0]->GetNumModes();
262  int order1 = m_base[1]->GetNumModes();
263 
264  Array<OneD, NekDouble> wsp(nquad2 * order0 * order1 +
265  nquad2 * nquad1 * nquad0);
266 
267  v_BwdTrans_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetBdata(),
268  m_base[2]->GetBdata(), inarray, outarray, wsp, true,
269  true, true);
270 }
271 
273  const Array<OneD, const NekDouble> &base0,
274  const Array<OneD, const NekDouble> &base1,
275  const Array<OneD, const NekDouble> &base2,
276  const Array<OneD, const NekDouble> &inarray,
278  bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
279 {
280  boost::ignore_unused(doCheckCollDir0, doCheckCollDir1, doCheckCollDir2);
281 
282  int nquad0 = m_base[0]->GetNumPoints();
283  int nquad1 = m_base[1]->GetNumPoints();
284  int nquad2 = m_base[2]->GetNumPoints();
285 
286  int order0 = m_base[0]->GetNumModes();
287  int order1 = m_base[1]->GetNumModes();
288  int order2 = m_base[2]->GetNumModes();
289 
290  Array<OneD, NekDouble> tmp = wsp;
291  Array<OneD, NekDouble> tmp1 = tmp + nquad2 * order0 * order1;
292 
293  int i, j, mode, mode1, cnt;
294 
295  // Perform summation over '2' direction
296  mode = mode1 = cnt = 0;
297  for (i = 0; i < order0; ++i)
298  {
299  for (j = 0; j < order1; ++j, ++cnt)
300  {
301  int ijmax = max(i, j);
302  Blas::Dgemv('N', nquad2, order2 - ijmax, 1.0,
303  base2.get() + mode * nquad2, nquad2,
304  inarray.get() + mode1, 1, 0.0, tmp.get() + cnt * nquad2,
305  1);
306  mode += order2 - ijmax;
307  mode1 += order2 - ijmax;
308  }
309  // increment mode in case order1!=order2
310  for (j = order1; j < order2 - i; ++j)
311  {
312  int ijmax = max(i, j);
313  mode += order2 - ijmax;
314  }
315  }
316 
317  // fix for modified basis by adding split of top singular
318  // vertex mode - currently (1+c)/2 x (1-b)/2 x (1-a)/2
319  // component is evaluated
321  {
322 
323  // Not sure why we could not use basis as 1.0
324  // top singular vertex - (1+c)/2 x (1+b)/2 x (1-a)/2 component
325  Blas::Daxpy(nquad2, inarray[1], base2.get() + nquad2, 1,
326  &tmp[0] + nquad2, 1);
327 
328  // top singular vertex - (1+c)/2 x (1-b)/2 x (1+a)/2 component
329  Blas::Daxpy(nquad2, inarray[1], base2.get() + nquad2, 1,
330  &tmp[0] + order1 * nquad2, 1);
331 
332  // top singular vertex - (1+c)/2 x (1+b)/2 x (1+a)/2 component
333  Blas::Daxpy(nquad2, inarray[1], base2.get() + nquad2, 1,
334  &tmp[0] + order1 * nquad2 + nquad2, 1);
335  }
336 
337  // Perform summation over '1' direction
338  mode = 0;
339  for (i = 0; i < order0; ++i)
340  {
341  Blas::Dgemm('N', 'T', nquad1, nquad2, order1, 1.0, base1.get(), nquad1,
342  tmp.get() + mode * nquad2, nquad2, 0.0,
343  tmp1.get() + i * nquad1 * nquad2, nquad1);
344  mode += order1;
345  }
346 
347  // Perform summation over '0' direction
348  Blas::Dgemm('N', 'T', nquad0, nquad1 * nquad2, order0, 1.0, base0.get(),
349  nquad0, tmp1.get(), nquad1 * nquad2, 0.0, outarray.get(),
350  nquad0);
351 }
352 
353 /** \brief Forward transform from physical quadrature space
354  stored in \a inarray and evaluate the expansion coefficients and
355  store in \a outarray
356 
357  Inputs:\n
358 
359  - \a inarray: array of physical quadrature points to be transformed
360 
361  Outputs:\n
362 
363  - \a outarray: updated array of expansion coefficients.
364 
365 */
367  Array<OneD, NekDouble> &outarray)
368 {
369  v_IProductWRTBase(inarray, outarray);
370 
371  // get Mass matrix inverse
372  StdMatrixKey imasskey(eInvMass, DetShapeType(), *this);
373  DNekMatSharedPtr imatsys = GetStdMatrix(imasskey);
374 
375  // copy inarray in case inarray == outarray
376  DNekVec in(m_ncoeffs, outarray);
377  DNekVec out(m_ncoeffs, outarray, eWrapper);
378 
379  out = (*imatsys) * in;
380 }
381 
382 //---------------------------------------
383 // Inner product functions
384 //---------------------------------------
385 
386 /** \brief Inner product of \a inarray over region with respect to the
387  expansion basis m_base[0]->GetBdata(),m_base[1]->GetBdata(),
388  m_base[2]->GetBdata() and return in \a outarray
389 
390  Wrapper call to StdPyrExp::IProductWRTBase
391 
392  Input:\n
393 
394  - \a inarray: array of function evaluated at the physical collocation points
395 
396  Output:\n
397 
398  - \a outarray: array of inner product with respect to each basis over region
399 
400 */
402  Array<OneD, NekDouble> &outarray)
403 {
404  if (m_base[0]->Collocation() && m_base[1]->Collocation() &&
405  m_base[2]->Collocation())
406  {
407  v_MultiplyByStdQuadratureMetric(inarray, outarray);
408  }
409  else
410  {
411  StdPyrExp::v_IProductWRTBase_SumFac(inarray, outarray);
412  }
413 }
414 
416  const Array<OneD, const NekDouble> &inarray,
417  Array<OneD, NekDouble> &outarray, bool multiplybyweights)
418 {
419 
420  int nquad1 = m_base[1]->GetNumPoints();
421  int nquad2 = m_base[2]->GetNumPoints();
422  int order0 = m_base[0]->GetNumModes();
423  int order1 = m_base[1]->GetNumModes();
424 
425  Array<OneD, NekDouble> wsp(order0 * nquad2 * (nquad1 + order1));
426 
427  if (multiplybyweights)
428  {
429  Array<OneD, NekDouble> tmp(inarray.size());
430 
431  v_MultiplyByStdQuadratureMetric(inarray, tmp);
432 
434  m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
435  tmp, outarray, wsp, true, true, true);
436  }
437  else
438  {
440  m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
441  inarray, outarray, wsp, true, true, true);
442  }
443 }
444 
446  const Array<OneD, const NekDouble> &base0,
447  const Array<OneD, const NekDouble> &base1,
448  const Array<OneD, const NekDouble> &base2,
449  const Array<OneD, const NekDouble> &inarray,
451  bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
452 {
453  boost::ignore_unused(doCheckCollDir0, doCheckCollDir1, doCheckCollDir2);
454 
455  int nquad0 = m_base[0]->GetNumPoints();
456  int nquad1 = m_base[1]->GetNumPoints();
457  int nquad2 = m_base[2]->GetNumPoints();
458 
459  int order0 = m_base[0]->GetNumModes();
460  int order1 = m_base[1]->GetNumModes();
461  int order2 = m_base[2]->GetNumModes();
462 
463  ASSERTL1(wsp.size() >= nquad1 * nquad2 * order0 + nquad2 * order0 * order1,
464  "Insufficient workspace size");
465 
466  Array<OneD, NekDouble> tmp1 = wsp;
467  Array<OneD, NekDouble> tmp2 = wsp + nquad1 * nquad2 * order0;
468 
469  int i, j, mode, mode1, cnt;
470 
471  // Inner product with respect to the '0' direction
472  Blas::Dgemm('T', 'N', nquad1 * nquad2, order0, nquad0, 1.0, inarray.get(),
473  nquad0, base0.get(), nquad0, 0.0, tmp1.get(), nquad1 * nquad2);
474 
475  // Inner product with respect to the '1' direction
476  for (mode = i = 0; i < order0; ++i)
477  {
478  Blas::Dgemm('T', 'N', nquad2, order1, nquad1, 1.0,
479  tmp1.get() + i * nquad1 * nquad2, nquad1, base1.get(),
480  nquad1, 0.0, tmp2.get() + mode * nquad2, nquad2);
481  mode += order1;
482  }
483 
484  // Inner product with respect to the '2' direction
485  mode = mode1 = cnt = 0;
486  for (i = 0; i < order0; ++i)
487  {
488  for (j = 0; j < order1; ++j, ++cnt)
489  {
490  int ijmax = max(i, j);
491 
492  Blas::Dgemv('T', nquad2, order2 - ijmax, 1.0,
493  base2.get() + mode * nquad2, nquad2,
494  tmp2.get() + cnt * nquad2, 1, 0.0,
495  outarray.get() + mode1, 1);
496  mode += order2 - ijmax;
497  mode1 += order2 - ijmax;
498  }
499 
500  // increment mode in case order1!=order2
501  for (j = order1; j < order2; ++j)
502  {
503  int ijmax = max(i, j);
504  mode += order2 - ijmax;
505  }
506  }
507 
508  // fix for modified basis for top singular vertex component
509  // Already have evaluated (1+c)/2 (1-b)/2 (1-a)/2
511  {
512  // add in (1+c)/2 (1+b)/2 (1-a)/2 component
513  outarray[1] +=
514  Blas::Ddot(nquad2, base2.get() + nquad2, 1, &tmp2[nquad2], 1);
515 
516  // add in (1+c)/2 (1-b)/2 (1+a)/2 component
517  outarray[1] += Blas::Ddot(nquad2, base2.get() + nquad2, 1,
518  &tmp2[nquad2 * order1], 1);
519 
520  // add in (1+c)/2 (1+b)/2 (1+a)/2 component
521  outarray[1] += Blas::Ddot(nquad2, base2.get() + nquad2, 1,
522  &tmp2[nquad2 * order1 + nquad2], 1);
523  }
524 }
525 
527  const int dir, const Array<OneD, const NekDouble> &inarray,
528  Array<OneD, NekDouble> &outarray)
529 {
530  StdPyrExp::v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
531 }
532 
533 /**
534  * @param inarray Function evaluated at physical collocation
535  * points.
536  * @param outarray Inner product with respect to each basis
537  * function over the element.
538  */
540  const int dir, const Array<OneD, const NekDouble> &inarray,
541  Array<OneD, NekDouble> &outarray)
542 {
543  int i;
544  int nquad0 = m_base[0]->GetNumPoints();
545  int nquad1 = m_base[1]->GetNumPoints();
546  int nquad2 = m_base[2]->GetNumPoints();
547  int nqtot = nquad0 * nquad1 * nquad2;
548 
549  Array<OneD, NekDouble> gfac0(nquad0);
550  Array<OneD, NekDouble> gfac1(nquad1);
551  Array<OneD, NekDouble> gfac2(nquad2);
552  Array<OneD, NekDouble> tmp0(nqtot);
553 
554  int order0 = m_base[0]->GetNumModes();
555  int order1 = m_base[1]->GetNumModes();
556 
557  Array<OneD, NekDouble> wsp(nquad1 * nquad2 * order0 +
558  nquad2 * order0 * order1);
559 
560  const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
561  const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
562  const Array<OneD, const NekDouble> &z2 = m_base[2]->GetZ();
563 
564  // set up geometric factor: (1+z0)/2
565  for (i = 0; i < nquad0; ++i)
566  {
567  gfac0[i] = 0.5 * (1 + z0[i]);
568  }
569 
570  // set up geometric factor: (1+z1)/2
571  for (i = 0; i < nquad1; ++i)
572  {
573  gfac1[i] = 0.5 * (1 + z1[i]);
574  }
575 
576  // Set up geometric factor: 2/(1-z2)
577  for (i = 0; i < nquad2; ++i)
578  {
579  gfac2[i] = 2.0 / (1 - z2[i]);
580  }
581 
582  // Derivative in first/second direction is always scaled as follows
583  const int nq01 = nquad0 * nquad1;
584  for (i = 0; i < nquad2; ++i)
585  {
586  Vmath::Smul(nq01, gfac2[i], &inarray[0] + i * nq01, 1,
587  &tmp0[0] + i * nq01, 1);
588  }
589 
591 
592  switch (dir)
593  {
594  case 0:
595  {
597  m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
598  m_base[2]->GetBdata(), tmp0, outarray, wsp, false, true, true);
599  break;
600  }
601  case 1:
602  {
604  m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
605  m_base[2]->GetBdata(), tmp0, outarray, wsp, true, false, true);
606  break;
607  }
608  case 2:
609  {
612 
613  // Scale eta_1 derivative by gfac0
614  for (i = 0; i < nquad1 * nquad2; ++i)
615  {
616  Vmath::Vmul(nquad0, tmp0.get() + i * nquad0, 1, gfac0.get(), 1,
617  tmp0.get() + i * nquad0, 1);
618  }
620  m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
621  m_base[2]->GetBdata(), tmp0, tmp3, wsp, false, true, true);
622 
623  // Scale eta_2 derivative by gfac1*gfac2
624  for (i = 0; i < nquad2; ++i)
625  {
626  Vmath::Smul(nq01, gfac2[i], &inarray[0] + i * nq01, 1,
627  &tmp0[0] + i * nq01, 1);
628  }
629  for (i = 0; i < nquad1 * nquad2; ++i)
630  {
631  Vmath::Smul(nquad0, gfac1[i % nquad1], &tmp0[0] + i * nquad0, 1,
632  &tmp0[0] + i * nquad0, 1);
633  }
634 
637  m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
638  m_base[2]->GetBdata(), tmp0, tmp4, wsp, true, false, true);
639 
640  v_MultiplyByStdQuadratureMetric(inarray, tmp0);
642  m_base[0]->GetBdata(), m_base[1]->GetBdata(),
643  m_base[2]->GetDbdata(), tmp0, outarray, wsp, true, true, false);
644 
645  Vmath::Vadd(m_ncoeffs, &tmp3[0], 1, &outarray[0], 1, &outarray[0],
646  1);
647  Vmath::Vadd(m_ncoeffs, &tmp4[0], 1, &outarray[0], 1, &outarray[0],
648  1);
649  break;
650  }
651  default:
652  {
653  ASSERTL1(false, "input dir is out of range");
654  break;
655  }
656  }
657 }
658 
659 //---------------------------------------
660 // Evaluation functions
661 //---------------------------------------
662 
665 {
666  NekDouble d2 = 1.0 - xi[2];
667  if (fabs(d2) < NekConstants::kNekZeroTol)
668  {
669  if (d2 >= 0.)
670  {
672  }
673  else
674  {
676  }
677  }
678  eta[2] = xi[2]; // eta_z = xi_z
679  eta[1] = 2.0 * (1.0 + xi[1]) / d2 - 1.0;
680  eta[0] = 2.0 * (1.0 + xi[0]) / d2 - 1.0;
681 }
682 
685 {
686  xi[0] = (1.0 + eta[0]) * (1.0 - eta[2]) * 0.5 - 1.0;
687  xi[1] = (1.0 + eta[1]) * (1.0 - eta[2]) * 0.5 - 1.0;
688  xi[2] = eta[2];
689 }
690 
694 {
695  Array<OneD, const NekDouble> etaBar_x = m_base[0]->GetZ();
696  Array<OneD, const NekDouble> eta_y = m_base[1]->GetZ();
697  Array<OneD, const NekDouble> eta_z = m_base[2]->GetZ();
698  int Qx = GetNumPoints(0);
699  int Qy = GetNumPoints(1);
700  int Qz = GetNumPoints(2);
701 
702  // Convert collapsed coordinates into cartesian coordinates: eta --> xi
703  for (int k = 0; k < Qz; ++k)
704  {
705  for (int j = 0; j < Qy; ++j)
706  {
707  for (int i = 0; i < Qx; ++i)
708  {
709  int s = i + Qx * (j + Qy * k);
710 
711  xi_z[s] = eta_z[k];
712  xi_y[s] = (1.0 + eta_y[j]) * (1.0 - eta_z[k]) / 2.0 - 1.0;
713  xi_x[s] = (1.0 + etaBar_x[i]) * (1.0 - eta_z[k]) / 2.0 - 1.0;
714  }
715  }
716  }
717 }
718 
720  const Array<OneD, const NekDouble> &inarray,
721  std::array<NekDouble, 3> &firstOrderDerivs)
722 {
723  // Collapse coordinates
724  Array<OneD, NekDouble> coll(3, 0.0);
725  LocCoordToLocCollapsed(coord, coll);
726 
727  // If near singularity do the old interpolation matrix method
728  if ((1 - coll[2]) < 1e-5)
729  {
730  int totPoints = GetTotPoints();
731  Array<OneD, NekDouble> EphysDeriv0(totPoints), EphysDeriv1(totPoints),
732  EphysDeriv2(totPoints);
733  PhysDeriv(inarray, EphysDeriv0, EphysDeriv1, EphysDeriv2);
734 
736  I[0] = GetBase()[0]->GetI(coll);
737  I[1] = GetBase()[1]->GetI(coll + 1);
738  I[2] = GetBase()[2]->GetI(coll + 2);
739 
740  firstOrderDerivs[0] = PhysEvaluate(I, EphysDeriv0);
741  firstOrderDerivs[1] = PhysEvaluate(I, EphysDeriv1);
742  firstOrderDerivs[2] = PhysEvaluate(I, EphysDeriv2);
743  return PhysEvaluate(I, inarray);
744  }
745 
746  std::array<NekDouble, 3> interDeriv;
747  NekDouble val = StdExpansion3D::BaryTensorDeriv(coll, inarray, interDeriv);
748 
749  NekDouble fac = 2.0 / (1.0 - coll[2]);
750 
751  firstOrderDerivs[0] = fac * interDeriv[0];
752  firstOrderDerivs[1] = fac * interDeriv[1];
753  firstOrderDerivs[2] = ((1.0 + coll[0]) / (1.0 - coll[2])) * interDeriv[0] +
754  ((1.0 + coll[1]) / (1.0 - coll[2])) * interDeriv[1] +
755  interDeriv[2];
756 
757  return val;
758 }
759 
760 void StdPyrExp::v_FillMode(const int mode, Array<OneD, NekDouble> &outarray)
761 {
763  tmp[mode] = 1.0;
764  v_BwdTrans(tmp, outarray);
765 }
766 
767 void StdPyrExp::v_GetTraceNumModes(const int fid, int &numModes0,
768  int &numModes1, Orientation faceOrient)
769 {
770  int nummodes[3] = {m_base[0]->GetNumModes(), m_base[1]->GetNumModes(),
771  m_base[2]->GetNumModes()};
772  switch (fid)
773  {
774  // quad
775  case 0:
776  {
777  numModes0 = nummodes[0];
778  numModes1 = nummodes[1];
779  }
780  break;
781  case 1:
782  case 3:
783  {
784  numModes0 = nummodes[0];
785  numModes1 = nummodes[2];
786  }
787  break;
788  case 2:
789  case 4:
790  {
791  numModes0 = nummodes[1];
792  numModes1 = nummodes[2];
793  }
794  break;
795  }
796 
797  if (faceOrient >= 9)
798  {
799  std::swap(numModes0, numModes1);
800  }
801 }
802 
804  const Array<OneD, const NekDouble> &coords, int mode)
805 {
806  Array<OneD, NekDouble> coll(3);
807  LocCoordToLocCollapsed(coords, coll);
808 
809  const int nm0 = m_base[0]->GetNumModes();
810  const int nm1 = m_base[1]->GetNumModes();
811  const int nm2 = m_base[2]->GetNumModes();
812 
813  int mode0 = 0, mode1 = 0, mode2 = 0, cnt = 0;
814 
815  bool found = false;
816  for (mode0 = 0; mode0 < nm0; ++mode0)
817  {
818  for (mode1 = 0; mode1 < nm1; ++mode1)
819  {
820  int maxpq = max(mode0, mode1);
821  for (mode2 = 0; mode2 < nm2 - maxpq; ++mode2, ++cnt)
822  {
823  if (cnt == mode)
824  {
825  found = true;
826  break;
827  }
828  }
829 
830  if (found)
831  {
832  break;
833  }
834  }
835 
836  if (found)
837  {
838  break;
839  }
840 
841  for (int j = nm1; j < nm2; ++j)
842  {
843  int ijmax = max(mode0, j);
844  mode2 += nm2 - ijmax;
845  }
846  }
847 
848  if (mode == 1 && m_base[0]->GetBasisType() == LibUtilities::eModified_A)
849  {
850  return StdExpansion::BaryEvaluateBasis<2>(coll[2], 1);
851  }
852  else
853  {
854  return StdExpansion::BaryEvaluateBasis<0>(coll[0], mode0) *
855  StdExpansion::BaryEvaluateBasis<1>(coll[1], mode1) *
856  StdExpansion::BaryEvaluateBasis<2>(coll[2], mode2);
857  }
858 }
859 
860 //---------------------------------------
861 // Helper functions
862 //---------------------------------------
863 
865 {
866  return 5;
867 }
868 
870 {
871  return 8;
872 }
873 
875 {
876  return 5;
877 }
878 
880 {
881  return LibUtilities::ePyramid;
882 }
883 
885 {
888  "BasisType is not a boundary interior form");
891  "BasisType is not a boundary interior form");
894  "BasisType is not a boundary interior form");
895 
896  int P = m_base[0]->GetNumModes();
897  int Q = m_base[1]->GetNumModes();
898  int R = m_base[2]->GetNumModes();
899 
901 }
902 
904 {
907  "BasisType is not a boundary interior form");
910  "BasisType is not a boundary interior form");
913  "BasisType is not a boundary interior form");
914 
915  int P = m_base[0]->GetNumModes() - 1;
916  int Q = m_base[1]->GetNumModes() - 1;
917  int R = m_base[2]->GetNumModes() - 1;
918 
919  return (P + 1) * (Q + 1) // 1 rect. face on base
920  + 2 * (R + 1) + P * (1 + 2 * R - P) // 2 tri. (P,R) faces
921  + 2 * (R + 1) + Q * (1 + 2 * R - Q); // 2 tri. (Q,R) faces
922 }
923 
924 int StdPyrExp::v_GetTraceNcoeffs(const int i) const
925 {
926  ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
927 
928  if (i == 0)
929  {
930  return GetBasisNumModes(0) * GetBasisNumModes(1);
931  }
932  else if (i == 1 || i == 3)
933  {
934  int P = GetBasisNumModes(0) - 1, Q = GetBasisNumModes(2) - 1;
935  return Q + 1 + (P * (1 + 2 * Q - P)) / 2;
936  }
937  else
938  {
939  int P = GetBasisNumModes(1) - 1, Q = GetBasisNumModes(2) - 1;
940  return Q + 1 + (P * (1 + 2 * Q - P)) / 2;
941  }
942 }
943 
944 int StdPyrExp::v_GetTraceIntNcoeffs(const int i) const
945 {
946  ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
947 
948  int P = m_base[0]->GetNumModes() - 1;
949  int Q = m_base[1]->GetNumModes() - 1;
950  int R = m_base[2]->GetNumModes() - 1;
951 
952  if (i == 0)
953  {
954  return (P - 1) * (Q - 1);
955  }
956  else if (i == 1 || i == 3)
957  {
958  return (P - 1) * (2 * (R - 1) - (P - 1) - 1) / 2;
959  }
960  else
961  {
962  return (Q - 1) * (2 * (R - 1) - (Q - 1) - 1) / 2;
963  }
964 }
965 
966 int StdPyrExp::v_GetTraceNumPoints(const int i) const
967 {
968  ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
969 
970  if (i == 0)
971  {
972  return m_base[0]->GetNumPoints() * m_base[1]->GetNumPoints();
973  }
974  else if (i == 1 || i == 3)
975  {
976  return m_base[0]->GetNumPoints() * m_base[2]->GetNumPoints();
977  }
978  else
979  {
980  return m_base[1]->GetNumPoints() * m_base[2]->GetNumPoints();
981  }
982 }
983 
984 int StdPyrExp::v_GetEdgeNcoeffs(const int i) const
985 {
986  ASSERTL2(i >= 0 && i <= 7, "edge id is out of range");
987 
988  if (i == 0 || i == 2)
989  {
990  return GetBasisNumModes(0);
991  }
992  else if (i == 1 || i == 3)
993  {
994  return GetBasisNumModes(1);
995  }
996  else
997  {
998  return GetBasisNumModes(2);
999  }
1000 }
1001 
1003  const int k) const
1004 {
1005  ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
1006  ASSERTL2(k >= 0 && k <= 1, "basis key id is out of range");
1007 
1008  switch (i)
1009  {
1010  case 0:
1011  {
1013  m_base[k]->GetNumPoints(),
1014  m_base[k]->GetNumModes());
1015  }
1016  case 1:
1017  case 3:
1018  {
1019  return EvaluateTriFaceBasisKey(k, m_base[2 * k]->GetBasisType(),
1020  m_base[2 * k]->GetNumPoints(),
1021  m_base[2 * k]->GetNumModes());
1022  }
1023  case 2:
1024  case 4:
1025  {
1026  return EvaluateTriFaceBasisKey(k, m_base[k + 1]->GetBasisType(),
1027  m_base[k + 1]->GetNumPoints(),
1028  m_base[k + 1]->GetNumModes());
1029  }
1030  }
1031 
1032  // Should never get here.
1034 }
1035 
1037  const std::vector<unsigned int> &nummodes, int &modes_offset)
1038 {
1040  nummodes[modes_offset], nummodes[modes_offset + 1],
1041  nummodes[modes_offset + 2]);
1042 
1043  modes_offset += 3;
1044  return nmodes;
1045 }
1046 
1047 int StdPyrExp::v_GetVertexMap(int vId, bool useCoeffPacking)
1048 {
1052  "Mapping not defined for this type of basis");
1053 
1054  int l = 0;
1055 
1056  if (useCoeffPacking == true) // follow packing of coefficients i.e q,r,p
1057  {
1058  switch (vId)
1059  {
1060  case 0:
1061  l = GetMode(0, 0, 0);
1062  break;
1063  case 1:
1064  l = GetMode(0, 0, 1);
1065  break;
1066  case 2:
1067  l = GetMode(0, 1, 0);
1068  break;
1069  case 3:
1070  l = GetMode(1, 0, 0);
1071  break;
1072  case 4:
1073  l = GetMode(1, 1, 0);
1074  break;
1075  default:
1076  ASSERTL0(false, "local vertex id must be between 0 and 4");
1077  }
1078  }
1079  else
1080  {
1081  switch (vId)
1082  {
1083  case 0:
1084  l = GetMode(0, 0, 0);
1085  break;
1086  case 1:
1087  l = GetMode(1, 0, 0);
1088  break;
1089  case 2:
1090  l = GetMode(1, 1, 0);
1091  break;
1092  case 3:
1093  l = GetMode(0, 1, 0);
1094  break;
1095  case 4:
1096  l = GetMode(0, 0, 1);
1097  break;
1098  default:
1099  ASSERTL0(false, "local vertex id must be between 0 and 4");
1100  }
1101  }
1102 
1103  return l;
1104 }
1105 
1107 {
1110  "BasisType is not a boundary interior form");
1113  "BasisType is not a boundary interior form");
1116  "BasisType is not a boundary interior form");
1117 
1118  int P = m_base[0]->GetNumModes() - 1, p;
1119  int Q = m_base[1]->GetNumModes() - 1, q;
1120  int R = m_base[2]->GetNumModes() - 1, r;
1121 
1122  int nIntCoeffs = m_ncoeffs - NumBndryCoeffs();
1123 
1124  if (outarray.size() != nIntCoeffs)
1125  {
1126  outarray = Array<OneD, unsigned int>(nIntCoeffs);
1127  }
1128 
1129  int idx = 0;
1130 
1131  // Loop over all interior modes.
1132  for (p = 2; p <= P; ++p)
1133  {
1134  for (q = 2; q <= Q; ++q)
1135  {
1136  int maxpq = max(p, q);
1137  for (r = 1; r <= R - maxpq; ++r)
1138  {
1139  outarray[idx++] = GetMode(p, q, r);
1140  }
1141  }
1142  }
1143 }
1144 
1146 {
1149  "BasisType is not a boundary interior form");
1152  "BasisType is not a boundary interior form");
1155  "BasisType is not a boundary interior form");
1156 
1157  int P = m_base[0]->GetNumModes() - 1, p;
1158  int Q = m_base[1]->GetNumModes() - 1, q;
1159  int R = m_base[2]->GetNumModes() - 1, r;
1160  int idx = 0;
1161 
1162  int nBnd = NumBndryCoeffs();
1163 
1164  if (maparray.size() != nBnd)
1165  {
1166  maparray = Array<OneD, unsigned int>(nBnd);
1167  }
1168 
1169  // Loop over all boundary modes (in ascending order).
1170  for (p = 0; p <= P; ++p)
1171  {
1172  // First two q-r planes are entirely boundary modes.
1173  if (p <= 1)
1174  {
1175  for (q = 0; q <= Q; ++q)
1176  {
1177  int maxpq = max(p, q);
1178  for (r = 0; r <= R - maxpq; ++r)
1179  {
1180  maparray[idx++] = GetMode(p, q, r);
1181  }
1182  }
1183  }
1184  else
1185  {
1186  // Remaining q-r planes contain boundary modes on the two
1187  // front and back sides and edges 0 2.
1188  for (q = 0; q <= Q; ++q)
1189  {
1190  if (q <= 1)
1191  {
1192  for (r = 0; r <= R - p; ++r)
1193  {
1194  maparray[idx++] = GetMode(p, q, r);
1195  }
1196  }
1197  else
1198  {
1199  maparray[idx++] = GetMode(p, q, 0);
1200  }
1201  }
1202  }
1203  }
1204 }
1205 
1206 void StdPyrExp::v_GetTraceCoeffMap(const unsigned int fid,
1207  Array<OneD, unsigned int> &maparray)
1208 {
1210  "Method only implemented if BasisType is identical"
1211  "in x and y directions");
1214  "Method only implemented for Modified_A BasisType"
1215  "(x and y direction) and ModifiedPyr_C BasisType (z "
1216  "direction)");
1217 
1218  int p, q, r, P = 0, Q = 0, idx = 0;
1219 
1220  int order0 = m_base[0]->GetNumModes();
1221  int order1 = m_base[1]->GetNumModes();
1222  int order2 = m_base[2]->GetNumModes();
1223 
1224  switch (fid)
1225  {
1226  case 0:
1227  P = order0;
1228  Q = order1;
1229  break;
1230  case 1:
1231  case 3:
1232  P = order0;
1233  Q = order2;
1234  break;
1235  case 2:
1236  case 4:
1237  P = order1;
1238  Q = order2;
1239  break;
1240  default:
1241  ASSERTL0(false, "fid must be between 0 and 4");
1242  }
1243 
1244  if (maparray.size() != P * Q)
1245  {
1246  maparray = Array<OneD, unsigned int>(P * Q);
1247  }
1248 
1249  // Set up ordering inside each 2D face. Also for triangular faces,
1250  // populate signarray.
1251  switch (fid)
1252  {
1253  case 0: // Bottom quad
1254 
1255  for (q = 0; q < Q; ++q)
1256  {
1257  for (p = 0; p < P; ++p)
1258  {
1259  maparray[q * P + p] = GetMode(p, q, 0);
1260  }
1261  }
1262  break;
1263 
1264  case 1: // Front triangle
1265  for (p = 0; p < P; ++p)
1266  {
1267  for (r = 0; r < Q - p; ++r)
1268  {
1269  maparray[idx++] = GetMode(p, 0, r);
1270  }
1271  }
1272  break;
1273 
1274  case 2: // Right triangle
1275  maparray[idx++] = GetMode(1, 0, 0);
1276  maparray[idx++] = GetMode(0, 0, 1);
1277  for (r = 1; r < Q - 1; ++r)
1278  {
1279  maparray[idx++] = GetMode(1, 0, r);
1280  }
1281 
1282  for (q = 1; q < P; ++q)
1283  {
1284  for (r = 0; r < Q - q; ++r)
1285  {
1286  maparray[idx++] = GetMode(1, q, r);
1287  }
1288  }
1289  break;
1290 
1291  case 3: // Rear triangle
1292  maparray[idx++] = GetMode(0, 1, 0);
1293  maparray[idx++] = GetMode(0, 0, 1);
1294  for (r = 1; r < Q - 1; ++r)
1295  {
1296  maparray[idx++] = GetMode(0, 1, r);
1297  }
1298 
1299  for (p = 1; p < P; ++p)
1300  {
1301  for (r = 0; r < Q - p; ++r)
1302  {
1303  maparray[idx++] = GetMode(p, 1, r);
1304  }
1305  }
1306  break;
1307 
1308  case 4: // Left triangle
1309  for (q = 0; q < P; ++q)
1310  {
1311  for (r = 0; r < Q - q; ++r)
1312  {
1313  maparray[idx++] = GetMode(0, q, r);
1314  }
1315  }
1316  break;
1317 
1318  default:
1319  ASSERTL0(false, "Face to element map unavailable.");
1320  }
1321 }
1322 
1323 void StdPyrExp::v_GetElmtTraceToTraceMap(const unsigned int fid,
1324  Array<OneD, unsigned int> &maparray,
1325  Array<OneD, int> &signarray,
1326  Orientation faceOrient, int P, int Q)
1327 {
1329  "Method only implemented if BasisType is identical"
1330  "in x and y directions");
1333  "Method only implemented for Modified_A BasisType"
1334  "(x and y direction) and ModifiedPyr_C BasisType (z "
1335  "direction)");
1336 
1337  int i, j, k, p, r, nFaceCoeffs;
1338  int nummodesA = 0, nummodesB = 0;
1339 
1340  int order0 = m_base[0]->GetNumModes();
1341  int order1 = m_base[1]->GetNumModes();
1342  int order2 = m_base[2]->GetNumModes();
1343 
1344  switch (fid)
1345  {
1346  case 0:
1347  nummodesA = order0;
1348  nummodesB = order1;
1349  break;
1350  case 1:
1351  case 3:
1352  nummodesA = order0;
1353  nummodesB = order2;
1354  break;
1355  case 2:
1356  case 4:
1357  nummodesA = order1;
1358  nummodesB = order2;
1359  break;
1360  default:
1361  ASSERTL0(false, "fid must be between 0 and 4");
1362  }
1363 
1364  if (P == -1)
1365  {
1366  P = nummodesA;
1367  Q = nummodesB;
1368  nFaceCoeffs = GetTraceNcoeffs(fid);
1369  }
1370  else if (fid > 0)
1371  {
1372  nFaceCoeffs = P * (2 * Q - P + 1) / 2;
1373  }
1374  else
1375  {
1376  nFaceCoeffs = P * Q;
1377  }
1378 
1379  // Allocate the map array and sign array; set sign array to ones (+)
1380  if (maparray.size() != nFaceCoeffs)
1381  {
1382  maparray = Array<OneD, unsigned int>(nFaceCoeffs);
1383  }
1384 
1385  if (signarray.size() != nFaceCoeffs)
1386  {
1387  signarray = Array<OneD, int>(nFaceCoeffs, 1);
1388  }
1389  else
1390  {
1391  fill(signarray.get(), signarray.get() + nFaceCoeffs, 1);
1392  }
1393 
1394  // triangular faces
1395  if (fid > 0)
1396  {
1397  // zero signmap and set maparray to zero if elemental
1398  // modes are not as large as face modesl
1399  int idx = 0;
1400  int cnt = 0;
1401  int minPA = min(nummodesA, P);
1402  int minQB = min(nummodesB, Q);
1403 
1404  for (j = 0; j < minPA; ++j)
1405  {
1406  // set maparray
1407  for (k = 0; k < minQB - j; ++k, ++cnt)
1408  {
1409  maparray[idx++] = cnt;
1410  }
1411 
1412  cnt += nummodesB - minQB;
1413 
1414  for (k = nummodesB - j; k < Q - j; ++k)
1415  {
1416  signarray[idx] = 0.0;
1417  maparray[idx++] = maparray[0];
1418  }
1419  }
1420 #if 0 // not required?
1421 
1422  for (j = minPA; j < nummodesA; ++j)
1423  {
1424  // set maparray
1425  for (k = 0; k < minQB-j; ++k, ++cnt)
1426  {
1427  maparray[idx++] = cnt;
1428  }
1429 
1430  cnt += nummodesB-minQB;
1431 
1432  for (k = nummodesB-j; k < Q-j; ++k)
1433  {
1434  signarray[idx] = 0.0;
1435  maparray[idx++] = maparray[0];
1436  }
1437  }
1438 #endif
1439  for (j = nummodesA; j < P; ++j)
1440  {
1441  for (k = 0; k < Q - j; ++k)
1442  {
1443  signarray[idx] = 0.0;
1444  maparray[idx++] = maparray[0];
1445  }
1446  }
1447 
1448  // Triangles only have one possible orientation (base
1449  // direction reversed); swap edge modes.
1450  if (faceOrient == eDir1BwdDir1_Dir2FwdDir2)
1451  {
1452  swap(maparray[0], maparray[Q]);
1453  for (i = 1; i < Q - 1; ++i)
1454  {
1455  swap(maparray[i + 1], maparray[Q + i]);
1456  }
1457 
1458  idx = 0;
1459  for (p = 0; p < P; ++p)
1460  {
1461  for (r = 0; r < Q - p; ++r, idx++)
1462  {
1463  if (p > 1)
1464  {
1465  signarray[idx] = p % 2 ? -1 : 1;
1466  }
1467  }
1468  }
1469  }
1470  }
1471  else
1472  {
1473 
1474  // Set up an array indexing for quads, since the ordering may need
1475  // to be transposed.
1476  Array<OneD, int> arrayindx(nFaceCoeffs, -1);
1477 
1478  for (i = 0; i < Q; i++)
1479  {
1480  for (j = 0; j < P; j++)
1481  {
1482  if (faceOrient < eDir1FwdDir2_Dir2FwdDir1)
1483  {
1484  arrayindx[i * P + j] = i * P + j;
1485  }
1486  else
1487  {
1488  arrayindx[i * P + j] = j * Q + i;
1489  }
1490  }
1491  }
1492 
1493  // zero signmap and set maparray to zero if elemental
1494  // modes are not as large as face modesl
1495  for (j = 0; j < P; ++j)
1496  {
1497  // set up default maparray
1498  for (k = 0; k < Q; k++)
1499  {
1500  maparray[arrayindx[j + k * P]] = j + k * nummodesA;
1501  }
1502 
1503  for (k = nummodesB; k < Q; ++k)
1504  {
1505  signarray[arrayindx[j + k * P]] = 0.0;
1506  maparray[arrayindx[j + k * P]] = maparray[0];
1507  }
1508  }
1509 
1510  for (j = nummodesA; j < P; ++j)
1511  {
1512  for (k = 0; k < Q; ++k)
1513  {
1514  signarray[arrayindx[j + k * P]] = 0.0;
1515  maparray[arrayindx[j + k * P]] = maparray[0];
1516  }
1517  }
1518 
1519  // The code below is exactly the same as that taken from
1520  // StdHexExp and reverses the 'b' and 'a' directions as
1521  // appropriate (1st and 2nd if statements respectively) in
1522  // quadrilateral faces.
1523  if (faceOrient == eDir1FwdDir1_Dir2BwdDir2 ||
1524  faceOrient == eDir1BwdDir1_Dir2BwdDir2 ||
1525  faceOrient == eDir1BwdDir2_Dir2FwdDir1 ||
1526  faceOrient == eDir1BwdDir2_Dir2BwdDir1)
1527  {
1528  if (faceOrient < eDir1FwdDir2_Dir2FwdDir1)
1529  {
1530  for (i = 3; i < Q; i += 2)
1531  {
1532  for (j = 0; j < P; j++)
1533  {
1534  signarray[arrayindx[i * P + j]] *= -1;
1535  }
1536  }
1537 
1538  for (i = 0; i < P; i++)
1539  {
1540  swap(maparray[i], maparray[i + P]);
1541  swap(signarray[i], signarray[i + P]);
1542  }
1543  }
1544  else
1545  {
1546  for (i = 0; i < Q; i++)
1547  {
1548  for (j = 3; j < P; j += 2)
1549  {
1550  signarray[arrayindx[i * P + j]] *= -1;
1551  }
1552  }
1553 
1554  for (i = 0; i < Q; i++)
1555  {
1556  swap(maparray[i], maparray[i + Q]);
1557  swap(signarray[i], signarray[i + Q]);
1558  }
1559  }
1560  }
1561 
1562  if (faceOrient == eDir1BwdDir1_Dir2FwdDir2 ||
1563  faceOrient == eDir1BwdDir1_Dir2BwdDir2 ||
1564  faceOrient == eDir1FwdDir2_Dir2BwdDir1 ||
1565  faceOrient == eDir1BwdDir2_Dir2BwdDir1)
1566  {
1567  if (faceOrient < eDir1FwdDir2_Dir2FwdDir1)
1568  {
1569  for (i = 0; i < Q; i++)
1570  {
1571  for (j = 3; j < P; j += 2)
1572  {
1573  signarray[arrayindx[i * P + j]] *= -1;
1574  }
1575  }
1576 
1577  for (i = 0; i < Q; i++)
1578  {
1579  swap(maparray[i * P], maparray[i * P + 1]);
1580  swap(signarray[i * P], signarray[i * P + 1]);
1581  }
1582  }
1583  else
1584  {
1585  for (i = 3; i < Q; i += 2)
1586  {
1587  for (j = 0; j < P; j++)
1588  {
1589  signarray[arrayindx[i * P + j]] *= -1;
1590  }
1591  }
1592 
1593  for (i = 0; i < P; i++)
1594  {
1595  swap(maparray[i * Q], maparray[i * Q + 1]);
1596  swap(signarray[i * Q], signarray[i * Q + 1]);
1597  }
1598  }
1599  }
1600  }
1601 }
1602 
1604  const int eid, Array<OneD, unsigned int> &maparray,
1605  Array<OneD, int> &signarray, const Orientation edgeOrient)
1606 {
1607  int i;
1608  bool signChange;
1609  const int P = m_base[0]->GetNumModes() - 1;
1610  const int Q = m_base[1]->GetNumModes() - 1;
1611  const int R = m_base[2]->GetNumModes() - 1;
1612  const int nEdgeIntCoeffs = v_GetEdgeNcoeffs(eid) - 2;
1613 
1614  if (maparray.size() != nEdgeIntCoeffs)
1615  {
1616  maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
1617  }
1618 
1619  if (signarray.size() != nEdgeIntCoeffs)
1620  {
1621  signarray = Array<OneD, int>(nEdgeIntCoeffs, 1);
1622  }
1623  else
1624  {
1625  fill(signarray.get(), signarray.get() + nEdgeIntCoeffs, 1);
1626  }
1627 
1628  // If edge is oriented backwards, change sign of modes which have
1629  // degree 2n+1, n >= 1.
1630  signChange = edgeOrient == eBackwards;
1631 
1632  switch (eid)
1633  {
1634  case 0:
1635  for (i = 2; i <= P; ++i)
1636  {
1637  maparray[i - 2] = GetMode(i, 0, 0);
1638  }
1639  break;
1640 
1641  case 1:
1642  for (i = 2; i <= Q; ++i)
1643  {
1644  maparray[i - 2] = GetMode(1, i, 0);
1645  }
1646  break;
1647  case 2:
1648  for (i = 2; i <= P; ++i)
1649  {
1650  maparray[i - 2] = GetMode(i, 1, 0);
1651  }
1652  break;
1653 
1654  case 3:
1655  for (i = 2; i <= Q; ++i)
1656  {
1657  maparray[i - 2] = GetMode(0, i, 0);
1658  }
1659  break;
1660  case 4:
1661  for (i = 2; i <= R; ++i)
1662  {
1663  maparray[i - 2] = GetMode(0, 0, i);
1664  }
1665  break;
1666 
1667  case 5:
1668  for (i = 1; i <= R - 1; ++i)
1669  {
1670  maparray[i - 1] = GetMode(1, 0, i);
1671  }
1672  break;
1673  case 6:
1674  for (i = 1; i <= R - 1; ++i)
1675  {
1676  maparray[i - 1] = GetMode(1, 1, i);
1677  }
1678  break;
1679 
1680  case 7:
1681  for (i = 1; i <= R - 1; ++i)
1682  {
1683  maparray[i - 1] = GetMode(0, 1, i);
1684  }
1685  break;
1686  default:
1687  ASSERTL0(false, "Edge not defined.");
1688  break;
1689  }
1690 
1691  if (signChange)
1692  {
1693  for (i = 1; i < nEdgeIntCoeffs; i += 2)
1694  {
1695  signarray[i] = -1;
1696  }
1697  }
1698 }
1699 
1701  const int fid, Array<OneD, unsigned int> &maparray,
1702  Array<OneD, int> &signarray, const Orientation faceOrient)
1703 {
1704  const int P = m_base[0]->GetNumModes() - 1;
1705  const int Q = m_base[1]->GetNumModes() - 1;
1706  const int R = m_base[2]->GetNumModes() - 1;
1707  const int nFaceIntCoeffs = v_GetTraceIntNcoeffs(fid);
1708  int p, q, r, idx = 0;
1709  int nummodesA = 0;
1710  int nummodesB = 0;
1711  int i, j;
1712 
1713  if (maparray.size() != nFaceIntCoeffs)
1714  {
1715  maparray = Array<OneD, unsigned int>(nFaceIntCoeffs);
1716  }
1717 
1718  if (signarray.size() != nFaceIntCoeffs)
1719  {
1720  signarray = Array<OneD, int>(nFaceIntCoeffs, 1);
1721  }
1722  else
1723  {
1724  fill(signarray.get(), signarray.get() + nFaceIntCoeffs, 1);
1725  }
1726 
1727  // Set up an array indexing for quad faces, since the ordering may
1728  // need to be transposed depending on orientation.
1729  Array<OneD, int> arrayindx(nFaceIntCoeffs);
1730  if (fid == 0)
1731  {
1732  nummodesA = P - 1;
1733  nummodesB = Q - 1;
1734 
1735  for (i = 0; i < nummodesB; i++)
1736  {
1737  for (j = 0; j < nummodesA; j++)
1738  {
1739  if (faceOrient < 9)
1740  {
1741  arrayindx[i * nummodesA + j] = i * nummodesA + j;
1742  }
1743  else
1744  {
1745  arrayindx[i * nummodesA + j] = j * nummodesB + i;
1746  }
1747  }
1748  }
1749  }
1750 
1751  switch (fid)
1752  {
1753  case 0: // Bottom quad
1754  for (q = 2; q <= Q; ++q)
1755  {
1756  for (p = 2; p <= P; ++p)
1757  {
1758  maparray[arrayindx[(q - 2) * nummodesA + (p - 2)]] =
1759  GetMode(p, q, 0);
1760  }
1761  }
1762  break;
1763  case 1: // Front triangle
1764  for (p = 2; p <= P; ++p)
1765  {
1766  for (r = 1; r <= R - p; ++r)
1767  {
1768  if ((int)faceOrient == 7)
1769  {
1770  signarray[idx] = p % 2 ? -1 : 1;
1771  }
1772  maparray[idx++] = GetMode(p, 0, r);
1773  }
1774  }
1775  break;
1776  case 2: // Right triangle
1777  for (q = 2; q <= Q; ++q)
1778  {
1779  for (r = 1; r <= R - q; ++r)
1780  {
1781  if ((int)faceOrient == 7)
1782  {
1783  signarray[idx] = q % 2 ? -1 : 1;
1784  }
1785  maparray[idx++] = GetMode(1, q, r);
1786  }
1787  }
1788  break;
1789 
1790  case 3: // Rear triangle
1791  for (p = 2; p <= P; ++p)
1792  {
1793  for (r = 1; r <= R - p; ++r)
1794  {
1795  if ((int)faceOrient == 7)
1796  {
1797  signarray[idx] = p % 2 ? -1 : 1;
1798  }
1799  maparray[idx++] = GetMode(p, 1, r);
1800  }
1801  }
1802  break;
1803 
1804  case 4: // Left triangle
1805  for (q = 2; q <= Q; ++q)
1806  {
1807  for (r = 1; r <= R - q; ++r)
1808  {
1809  if ((int)faceOrient == 7)
1810  {
1811  signarray[idx] = q % 2 ? -1 : 1;
1812  }
1813  maparray[idx++] = GetMode(0, q, r);
1814  }
1815  }
1816  break;
1817  default:
1818  ASSERTL0(false, "Face interior map not available.");
1819  }
1820 
1821  // Triangular faces are processed in the above switch loop; for
1822  // remaining quad faces, set up orientation if necessary.
1823  if (fid > 0)
1824  {
1825  return;
1826  }
1827 
1828  if (faceOrient == 6 || faceOrient == 8 || faceOrient == 11 ||
1829  faceOrient == 12)
1830  {
1831  if (faceOrient < 9)
1832  {
1833  for (i = 1; i < nummodesB; i += 2)
1834  {
1835  for (j = 0; j < nummodesA; j++)
1836  {
1837  signarray[arrayindx[i * nummodesA + j]] *= -1;
1838  }
1839  }
1840  }
1841  else
1842  {
1843  for (i = 0; i < nummodesB; i++)
1844  {
1845  for (j = 1; j < nummodesA; j += 2)
1846  {
1847  signarray[arrayindx[i * nummodesA + j]] *= -1;
1848  }
1849  }
1850  }
1851  }
1852 
1853  if (faceOrient == 7 || faceOrient == 8 || faceOrient == 10 ||
1854  faceOrient == 12)
1855  {
1856  if (faceOrient < 9)
1857  {
1858  for (i = 0; i < nummodesB; i++)
1859  {
1860  for (j = 1; j < nummodesA; j += 2)
1861  {
1862  signarray[arrayindx[i * nummodesA + j]] *= -1;
1863  }
1864  }
1865  }
1866  else
1867  {
1868  for (i = 1; i < nummodesB; i += 2)
1869  {
1870  for (j = 0; j < nummodesA; j++)
1871  {
1872  signarray[arrayindx[i * nummodesA + j]] *= -1;
1873  }
1874  }
1875  }
1876  }
1877 }
1878 
1879 //---------------------------------------
1880 // Wrapper functions
1881 //---------------------------------------
1882 
1884 {
1885  return CreateGeneralMatrix(mkey);
1886 }
1887 
1889 {
1890  return v_GenMatrix(mkey);
1891 }
1892 
1893 /**
1894  * @brief Compute the mode number in the expansion for a
1895  * particular tensorial combination.
1896  *
1897  * Modes are numbered with the r index travelling fastest,
1898  * followed by q and then p, and each q-r plane is of size
1899  *
1900  * (R+1-p)*(Q+1) - l(l+1)/2 where l = max(0,Q-p)
1901  *
1902  * For example, when P=2, Q=3 and R=4 the indexing inside each
1903  * q-r plane (with r increasing upwards and q to the right)
1904  * is:
1905  *
1906  * p = 0: p = 1: p = 2:
1907  * ----------------------------------
1908  * 4
1909  * 3 8 17 21
1910  * 2 7 11 16 20 24 29 32 35
1911  * 1 6 10 13 15 19 23 26 28 31 34 37
1912  * 0 5 9 12 14 18 22 25 27 30 33 36
1913  *
1914  * Note that in this element, we must have that \f$ P,Q \leq
1915  * R\f$.
1916  */
1917 int StdPyrExp::GetMode(const int I, const int J, const int K)
1918 {
1919  const int Q = m_base[1]->GetNumModes() - 1;
1920  const int R = m_base[2]->GetNumModes() - 1;
1921 
1922  int i, l;
1923  int cnt = 0;
1924 
1925  // Traverse to q-r plane number I
1926  for (i = 0; i < I; ++i)
1927  {
1928  // Size of triangle part
1929  l = max(0, Q - i);
1930 
1931  // Size of rectangle part
1932  cnt += (R + 1 - i) * (Q + 1) - l * (l + 1) / 2;
1933  }
1934 
1935  // Traverse to q column J (Pretend this is a face of width J)
1936  l = max(0, J - 1 - I);
1937  cnt += (R + 1 - I) * J - l * (l + 1) / 2;
1938 
1939  // Traverse up stacks to K
1940  cnt += K;
1941 
1942  return cnt;
1943 }
1944 
1946  const Array<OneD, const NekDouble> &inarray,
1947  Array<OneD, NekDouble> &outarray)
1948 {
1949  int i, j;
1950 
1951  int nquad0 = m_base[0]->GetNumPoints();
1952  int nquad1 = m_base[1]->GetNumPoints();
1953  int nquad2 = m_base[2]->GetNumPoints();
1954 
1955  const Array<OneD, const NekDouble> &w0 = m_base[0]->GetW();
1956  const Array<OneD, const NekDouble> &w1 = m_base[1]->GetW();
1957  const Array<OneD, const NekDouble> &w2 = m_base[2]->GetW();
1958 
1959  const Array<OneD, const NekDouble> &z2 = m_base[2]->GetZ();
1960 
1961  // Multiply by integration constants in x-direction
1962  for (i = 0; i < nquad1 * nquad2; ++i)
1963  {
1964  Vmath::Vmul(nquad0, inarray.get() + i * nquad0, 1, w0.get(), 1,
1965  outarray.get() + i * nquad0, 1);
1966  }
1967 
1968  // Multiply by integration constants in y-direction
1969  for (j = 0; j < nquad2; ++j)
1970  {
1971  for (i = 0; i < nquad1; ++i)
1972  {
1973  Blas::Dscal(nquad0, w1[i],
1974  &outarray[0] + i * nquad0 + j * nquad0 * nquad1, 1);
1975  }
1976  }
1977 
1978  // Multiply by integration constants in z-direction; need to
1979  // incorporate factor [(1-eta_3)/2]^2 into weights, but only if
1980  // using GLL quadrature points.
1981  switch (m_base[2]->GetPointsType())
1982  {
1983  // (2,0) Jacobi inner product.
1984  case LibUtilities::eGaussRadauMAlpha2Beta0:
1985  for (i = 0; i < nquad2; ++i)
1986  {
1987  Blas::Dscal(nquad0 * nquad1, 0.25 * w2[i],
1988  &outarray[0] + i * nquad0 * nquad1, 1);
1989  }
1990  break;
1991 
1992  default:
1993  for (i = 0; i < nquad2; ++i)
1994  {
1995  Blas::Dscal(nquad0 * nquad1,
1996  0.25 * (1 - z2[i]) * (1 - z2[i]) * w2[i],
1997  &outarray[0] + i * nquad0 * nquad1, 1);
1998  }
1999  break;
2000  }
2001 }
2002 
2004  const StdMatrixKey &mkey)
2005 {
2006  // Generate an orthonogal expansion
2007  int qa = m_base[0]->GetNumPoints();
2008  int qb = m_base[1]->GetNumPoints();
2009  int qc = m_base[2]->GetNumPoints();
2010  int nmodes_a = m_base[0]->GetNumModes();
2011  int nmodes_b = m_base[1]->GetNumModes();
2012  int nmodes_c = m_base[2]->GetNumModes();
2013  // Declare orthogonal basis.
2017 
2021  StdPyrExp OrthoExp(Ba, Bb, Bc);
2022 
2023  Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
2024  int i, j, k, cnt = 0;
2025 
2026  // project onto modal space.
2027  OrthoExp.FwdTrans(array, orthocoeffs);
2028 
2030  {
2031  // Rodrigo's power kernel
2033  NekDouble SvvDiffCoeff =
2036 
2037  for (int i = 0; i < nmodes_a; ++i)
2038  {
2039  for (int j = 0; j < nmodes_b; ++j)
2040  {
2041  int maxij = max(i, j);
2042  NekDouble fac1 = std::max(
2043  pow((1.0 * i) / (nmodes_a - 1), cutoff * nmodes_a),
2044  pow((1.0 * j) / (nmodes_b - 1), cutoff * nmodes_b));
2045 
2046  for (int k = 0; k < nmodes_c - maxij; ++k)
2047  {
2048  NekDouble fac =
2049  std::max(fac1, pow((1.0 * k) / (nmodes_c - 1),
2050  cutoff * nmodes_c));
2051 
2052  orthocoeffs[cnt] *= SvvDiffCoeff * fac;
2053  cnt++;
2054  }
2055  }
2056  }
2057  }
2058  else if (mkey.ConstFactorExists(
2059  eFactorSVVDGKerDiffCoeff)) // Rodrigo/Mansoor's DG Kernel
2060  {
2061  NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDGKerDiffCoeff) *
2063 
2064  int max_abc = max(nmodes_a - kSVVDGFiltermodesmin,
2065  nmodes_b - kSVVDGFiltermodesmin);
2066  max_abc = max(max_abc, nmodes_c - kSVVDGFiltermodesmin);
2067  // clamp max_abc
2068  max_abc = max(max_abc, 0);
2069  max_abc = min(max_abc, kSVVDGFiltermodesmax - kSVVDGFiltermodesmin);
2070 
2071  for (int i = 0; i < nmodes_a; ++i)
2072  {
2073  for (int j = 0; j < nmodes_b; ++j)
2074  {
2075  int maxij = max(i, j);
2076 
2077  for (int k = 0; k < nmodes_c - maxij; ++k)
2078  {
2079  int maxijk = max(maxij, k);
2080  maxijk = min(maxijk, kSVVDGFiltermodesmax - 1);
2081 
2082  orthocoeffs[cnt] *=
2083  SvvDiffCoeff * kSVVDGFilter[max_abc][maxijk];
2084  cnt++;
2085  }
2086  }
2087  }
2088  }
2089  else
2090  {
2091  // SVV filter paramaters (how much added diffusion relative
2092  // to physical one and fraction of modes from which you
2093  // start applying this added diffusion)
2094  //
2095  NekDouble SvvDiffCoeff =
2097  NekDouble SVVCutOff =
2099 
2100  // Defining the cut of mode
2101  int cutoff_a = (int)(SVVCutOff * nmodes_a);
2102  int cutoff_b = (int)(SVVCutOff * nmodes_b);
2103  int cutoff_c = (int)(SVVCutOff * nmodes_c);
2104  // To avoid the fac[j] from blowing up
2105  NekDouble epsilon = 1;
2106 
2107  int nmodes = min(min(nmodes_a, nmodes_b), nmodes_c);
2108  NekDouble cutoff = min(min(cutoff_a, cutoff_b), cutoff_c);
2109 
2110  for (i = 0; i < nmodes_a; ++i) // P
2111  {
2112  for (j = 0; j < nmodes_b; ++j) // Q
2113  {
2114  int maxij = max(i, j);
2115  for (k = 0; k < nmodes_c - maxij; ++k) // R
2116  {
2117  if (j + k >= cutoff || i + k >= cutoff)
2118  {
2119  orthocoeffs[cnt] *=
2120  (SvvDiffCoeff *
2121  exp(-(i + k - nmodes) * (i + k - nmodes) /
2122  ((NekDouble)((i + k - cutoff + epsilon) *
2123  (i + k - cutoff + epsilon)))) *
2124  exp(-(j - nmodes) * (j - nmodes) /
2125  ((NekDouble)((j - cutoff + epsilon) *
2126  (j - cutoff + epsilon)))));
2127  }
2128  else
2129  {
2130  orthocoeffs[cnt] *= 0.0;
2131  }
2132  cnt++;
2133  }
2134  }
2135  }
2136  }
2137 
2138  // backward transform to physical space
2139  OrthoExp.BwdTrans(orthocoeffs, array);
2140 }
2141 
2143  const Array<OneD, const NekDouble> &inarray,
2144  Array<OneD, NekDouble> &outarray)
2145 {
2146  int nquad0 = m_base[0]->GetNumPoints();
2147  int nquad1 = m_base[1]->GetNumPoints();
2148  int nquad2 = m_base[2]->GetNumPoints();
2149  int nqtot = nquad0 * nquad1 * nquad2;
2150  int nmodes0 = m_base[0]->GetNumModes();
2151  int nmodes1 = m_base[1]->GetNumModes();
2152  int nmodes2 = m_base[2]->GetNumModes();
2153  int numMax = nmodes0;
2154 
2156  Array<OneD, NekDouble> coeff_tmp1(m_ncoeffs, 0.0);
2157  Array<OneD, NekDouble> phys_tmp(nqtot, 0.0);
2158  Array<OneD, NekDouble> tmp, tmp2, tmp3, tmp4;
2159 
2160  const LibUtilities::PointsKey Pkey0 = m_base[0]->GetPointsKey();
2161  const LibUtilities::PointsKey Pkey1 = m_base[1]->GetPointsKey();
2162  const LibUtilities::PointsKey Pkey2 = m_base[2]->GetPointsKey();
2163 
2164  LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A, nmodes0, Pkey0);
2165  LibUtilities::BasisKey bortho1(LibUtilities::eOrtho_A, nmodes1, Pkey1);
2166  LibUtilities::BasisKey bortho2(LibUtilities::eOrthoPyr_C, nmodes2, Pkey2);
2167 
2168  int cnt = 0;
2169  int u = 0;
2170  int i = 0;
2171  StdRegions::StdPyrExpSharedPtr OrthoPyrExp;
2172 
2174  bortho0, bortho1, bortho2);
2175 
2176  BwdTrans(inarray, phys_tmp);
2177  OrthoPyrExp->FwdTrans(phys_tmp, coeff);
2178 
2179  // filtering
2180  for (u = 0; u < numMin; ++u)
2181  {
2182  for (i = 0; i < numMin; ++i)
2183  {
2184 
2185  int maxui = max(u, i);
2186  Vmath::Vcopy(numMin - maxui, tmp = coeff + cnt, 1,
2187  tmp2 = coeff_tmp1 + cnt, 1);
2188  cnt += nmodes2 - maxui;
2189  }
2190 
2191  for (i = numMin; i < nmodes1; ++i)
2192  {
2193  int maxui = max(u, i);
2194  cnt += numMax - maxui;
2195  }
2196  }
2197 
2198  OrthoPyrExp->BwdTrans(coeff_tmp1, phys_tmp);
2199  StdPyrExp::FwdTrans(phys_tmp, outarray);
2200 }
2201 
2202 } // namespace StdRegions
2203 } // 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
Describes the specification for a Basis.
Definition: Basis.h:50
BasisType GetBasisType() const
Return type of expansion basis.
Definition: Basis.h:141
int GetNumModes() const
Returns the order of the basis.
Definition: Basis.h:83
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 IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
NekDouble BaryTensorDeriv(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d1, Array< OneD, NekDouble > &outarray_d2, Array< OneD, NekDouble > &outarray_d3)
Calculate the 3D derivative in the local tensor/collapsed coordinate at the physical points.
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
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140
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.
const Array< OneD, const LibUtilities::BasisSharedPtr > & GetBase() const
This function gets the shared point to basis.
Definition: StdExpansion.h:106
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
NekDouble PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
This function evaluates the expansion at a single (arbitrary) point of the domain.
Definition: StdExpansion.h:918
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:373
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
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.
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:224
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:175
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:848
Array< OneD, LibUtilities::BasisSharedPtr > m_base
NekDouble GetConstFactor(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:126
bool ConstFactorExists(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:135
void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdPyrExp.cpp:539
DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override
Definition: StdPyrExp.cpp:1883
int v_NumBndryCoeffs() const override
Definition: StdPyrExp.cpp:884
NekDouble v_PhysEvaluateBasis(const Array< OneD, const NekDouble > &coords, int mode) final
Definition: StdPyrExp.cpp:803
void v_ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdPyrExp.cpp:2142
int v_NumDGBndryCoeffs() const override
Definition: StdPyrExp.cpp:903
void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Forward transform from physical quadrature space stored in inarray and evaluate the expansion coeffic...
Definition: StdPyrExp.cpp:366
void v_GetElmtTraceToTraceMap(const unsigned int fid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation faceOrient, int P, int Q) override
Definition: StdPyrExp.cpp:1323
int v_GetTraceNcoeffs(const int i) const override
Definition: StdPyrExp.cpp:924
void v_LocCollapsedToLocCoord(const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi) override
Definition: StdPyrExp.cpp:683
void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey) override
Definition: StdPyrExp.cpp:2003
int v_GetEdgeNcoeffs(const int i) const override
Definition: StdPyrExp.cpp:984
void v_BwdTrans_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2) override
Definition: StdPyrExp.cpp:272
int v_GetNtraces() const override
Definition: StdPyrExp.cpp:874
void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdPyrExp.cpp:526
void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray) override
Definition: StdPyrExp.cpp:760
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
Definition: StdPyrExp.cpp:415
int v_GetNedges() const override
Definition: StdPyrExp.cpp:869
int v_CalcNumberOfCoefficients(const std::vector< unsigned int > &nummodes, int &modes_offset) override
Definition: StdPyrExp.cpp:1036
int GetMode(int I, int J, int K)
Compute the mode number in the expansion for a particular tensorial combination.
Definition: StdPyrExp.cpp:1917
int v_GetVertexMap(int localVertexId, bool useCoeffPacking=false) override
Definition: StdPyrExp.cpp:1047
int v_GetTraceNumPoints(const int i) const override
Definition: StdPyrExp.cpp:966
void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Backward transformation is evaluated at the quadrature points.
Definition: StdPyrExp.cpp:236
void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray) override
Definition: StdPyrExp.cpp:1145
void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Inner product of inarray over region with respect to the expansion basis m_base[0]->GetBdata(),...
Definition: StdPyrExp.cpp:401
void v_GetInteriorMap(Array< OneD, unsigned int > &outarray) override
Definition: StdPyrExp.cpp:1106
void v_GetEdgeInteriorToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2) override
Definition: StdPyrExp.cpp:1603
void v_LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta) override
Definition: StdPyrExp.cpp:663
DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey) override
Definition: StdPyrExp.cpp:1888
int v_GetNverts() const override
Definition: StdPyrExp.cpp:864
void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdPyrExp.cpp:1945
void v_GetTraceInteriorToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2) override
Definition: StdPyrExp.cpp:1700
void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdPyrExp.cpp:255
void v_IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2) override
Definition: StdPyrExp.cpp:445
void v_StdPhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2) override
Definition: StdPyrExp.cpp:196
const LibUtilities::BasisKey v_GetTraceBasisKey(const int i, const int k) const override
Definition: StdPyrExp.cpp:1002
LibUtilities::ShapeType v_DetShapeType() const override
Definition: StdPyrExp.cpp:879
void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2) override
Calculate the derivative of the physical points.
Definition: StdPyrExp.cpp:91
void v_GetTraceNumModes(const int fid, int &numModes0, int &numModes1, Orientation faceOrient=eDir1FwdDir1_Dir2FwdDir2) override
Definition: StdPyrExp.cpp:767
NekDouble v_PhysEvaluate(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
Definition: StdPyrExp.cpp:719
void v_GetTraceCoeffMap(const unsigned int fid, Array< OneD, unsigned int > &maparray) override
Definition: StdPyrExp.cpp:1206
void v_GetCoords(Array< OneD, NekDouble > &xi_x, Array< OneD, NekDouble > &xi_y, Array< OneD, NekDouble > &xi_z) override
Definition: StdPyrExp.cpp:691
int v_GetTraceIntNcoeffs(const int i) const override
Definition: StdPyrExp.cpp:944
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 Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition: Blas.hpp:168
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
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:154
int getNumberOfCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:237
int getNumberOfBndCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:264
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:44
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:58
@ eModifiedPyr_C
Principle Modified Functions.
Definition: BasisType.h:55
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:50
@ eOrthoPyr_C
Principle Orthogonal Functions .
Definition: BasisType.h:53
static const NekDouble kNekZeroTol
LibUtilities::BasisKey EvaluateTriFaceBasisKey(const int facedir, const LibUtilities::BasisType faceDirBasisType, const int numpoints, const int nummodes)
const int kSVVDGFiltermodesmin
Definition: StdRegions.hpp:469
std::shared_ptr< StdPyrExp > StdPyrExpSharedPtr
Definition: StdPyrExp.h:258
const int kSVVDGFiltermodesmax
Definition: StdRegions.hpp:470
LibUtilities::BasisKey EvaluateQuadFaceBasisKey(const int facedir, const LibUtilities::BasisType faceDirBasisType, const int numpoints, const int nummodes)
const NekDouble kSVVDGFilter[9][11]
Definition: StdRegions.hpp:472
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 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:359
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 Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255