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