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