Nektar++
Loading...
Searching...
No Matches
StdTetExp.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: StdTetExp.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: Header field for tetrahedral routines built upon
32// StdExpansion3D
33//
34///////////////////////////////////////////////////////////////////////////////
35
37
38using namespace std;
39
40namespace Nektar::StdRegions
41{
43 const LibUtilities::BasisKey &Bb,
44 const LibUtilities::BasisKey &Bc)
45 : StdExpansion(LibUtilities::StdTetData::getNumberOfCoefficients(
46 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
47 3, Ba, Bb, Bc),
48 StdExpansion3D(LibUtilities::StdTetData::getNumberOfCoefficients(
49 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
50 Ba, Bb, Bc)
51{
52 ASSERTL0(Ba.GetNumModes() <= Bb.GetNumModes(),
53 "order in 'a' direction is higher than order "
54 "in 'b' direction");
55 ASSERTL0(Ba.GetNumModes() <= Bc.GetNumModes(),
56 "order in 'a' direction is higher than order "
57 "in 'c' direction");
58 ASSERTL0(Bb.GetNumModes() <= Bc.GetNumModes(),
59 "order in 'b' direction is higher than order "
60 "in 'c' direction");
61}
62
63//----------------------------
64// Differentiation Methods
65//----------------------------
66
67/**
68 * \brief Calculate the derivative of the physical points
69 *
70 * The derivative is evaluated at the nodal physical points.
71 * Derivatives with respect to the local Cartesian coordinates
72 *
73 * \f$\begin{Bmatrix} \frac {\partial} {\partial \xi_1} \\ \frac
74 * {\partial} {\partial \xi_2} \\ \frac {\partial} {\partial \xi_3}
75 * \end{Bmatrix} = \begin{Bmatrix} \frac 4 {(1-\eta_2)(1-\eta_3)}
76 * \frac \partial {\partial \eta_1} \ \ \frac {2(1+\eta_1)}
77 * {(1-\eta_2)(1-\eta_3)} \frac \partial {\partial \eta_1} + \frac 2
78 * {1-\eta_3} \frac \partial {\partial \eta_3} \\ \frac {2(1 +
79 * \eta_1)} {2(1 - \eta_2)(1-\eta_3)} \frac \partial {\partial \eta_1}
80 * + \frac {1 + \eta_2} {1 - \eta_3} \frac \partial {\partial \eta_2}
81 * + \frac \partial {\partial \eta_3} \end{Bmatrix}\f$
82 **/
84 Array<OneD, NekDouble> &out_dxi0,
85 Array<OneD, NekDouble> &out_dxi1,
86 Array<OneD, NekDouble> &out_dxi2)
87{
88 int Q0 = m_base[0]->GetNumPoints();
89 int Q1 = m_base[1]->GetNumPoints();
90 int Q2 = m_base[2]->GetNumPoints();
91 int Qtot = Q0 * Q1 * Q2;
92
93 // Compute the physical derivative
94 Array<OneD, NekDouble> out_dEta0(3 * Qtot, 0.0);
95 Array<OneD, NekDouble> out_dEta1 = out_dEta0 + Qtot;
96 Array<OneD, NekDouble> out_dEta2 = out_dEta1 + Qtot;
97
98 bool Do_2 = (out_dxi2.size() > 0) ? true : false;
99 bool Do_1 = (out_dxi1.size() > 0) ? true : false;
100
101 if (Do_2) // Need all local derivatives
102 {
103 PhysTensorDeriv(inarray, out_dEta0, out_dEta1, out_dEta2);
104 }
105 else if (Do_1) // Need 0 and 1 derivatives
106 {
107 PhysTensorDeriv(inarray, out_dEta0, out_dEta1, NullNekDouble1DArray);
108 }
109 else // Only need Eta0 derivaitve
110 {
111 PhysTensorDeriv(inarray, out_dEta0, NullNekDouble1DArray,
113 }
114
115 Array<OneD, const NekDouble> eta_0, eta_1, eta_2;
116 eta_0 = m_base[0]->GetZ();
117 eta_1 = m_base[1]->GetZ();
118 eta_2 = m_base[2]->GetZ();
119
120 // calculate 2.0/((1-eta_1)(1-eta_2)) Out_dEta0
121
122 NekDouble *dEta0 = &out_dEta0[0];
123 NekDouble fac;
124 for (int k = 0; k < Q2; ++k)
125 {
126 for (int j = 0; j < Q1; ++j, dEta0 += Q0)
127 {
128 Vmath::Smul(Q0, 2.0 / (1.0 - eta_1[j]), dEta0, 1, dEta0, 1);
129 }
130 fac = 1.0 / (1.0 - eta_2[k]);
131 Vmath::Smul(Q0 * Q1, fac, &out_dEta0[0] + k * Q0 * Q1, 1,
132 &out_dEta0[0] + k * Q0 * Q1, 1);
133 }
134
135 if (out_dxi0.size() > 0)
136 {
137 // out_dxi0 = 4.0/((1-eta_1)(1-eta_2)) Out_dEta0
138 Vmath::Smul(Qtot, 2.0, out_dEta0, 1, out_dxi0, 1);
139 }
140
141 if (Do_1 || Do_2)
142 {
143 Array<OneD, NekDouble> Fac0(Q0);
144 Vmath::Sadd(Q0, 1.0, eta_0, 1, Fac0, 1);
145
146 // calculate 2.0*(1+eta_0)/((1-eta_1)(1-eta_2)) Out_dEta0
147 for (int k = 0; k < Q1 * Q2; ++k)
148 {
149 Vmath::Vmul(Q0, &Fac0[0], 1, &out_dEta0[0] + k * Q0, 1,
150 &out_dEta0[0] + k * Q0, 1);
151 }
152 // calculate 2/(1.0-eta_2) out_dEta1
153 for (int k = 0; k < Q2; ++k)
154 {
155 Vmath::Smul(Q0 * Q1, 2.0 / (1.0 - eta_2[k]),
156 &out_dEta1[0] + k * Q0 * Q1, 1,
157 &out_dEta1[0] + k * Q0 * Q1, 1);
158 }
159
160 if (Do_1)
161 {
162 // calculate out_dxi1 = 2.0(1+eta_0)/((1-eta_1)(1-eta_2)) Out_dEta0
163 // + 2/(1.0-eta_2) out_dEta1
164 Vmath::Vadd(Qtot, out_dEta0, 1, out_dEta1, 1, out_dxi1, 1);
165 }
166
167 if (Do_2)
168 {
169 // calculate (1 + eta_1)/(1 -eta_2)*out_dEta1
170 NekDouble *dEta1 = &out_dEta1[0];
171 for (int k = 0; k < Q2; ++k)
172 {
173 for (int j = 0; j < Q1; ++j, dEta1 += Q0)
174 {
175 Vmath::Smul(Q0, (1.0 + eta_1[j]) / 2.0, dEta1, 1, dEta1, 1);
176 }
177 }
178
179 // calculate out_dxi2 =
180 // 2.0(1+eta_0)/((1-eta_1)(1-eta_2)) Out_dEta0 +
181 // (1 + eta_1)/(1 -eta_2)*out_dEta1 + out_dEta2
182 Vmath::Vadd(Qtot, out_dEta0, 1, out_dEta1, 1, out_dxi2, 1);
183 Vmath::Vadd(Qtot, out_dEta2, 1, out_dxi2, 1, out_dxi2, 1);
184 }
185 }
186}
187
188/**
189 * @param dir Direction in which to compute derivative.
190 * Valid values are 0, 1, 2.
191 * @param inarray Input array.
192 * @param outarray Output array.
193 */
194void StdTetExp::v_PhysDeriv(const int dir,
195 const Array<OneD, const NekDouble> &inarray,
196 Array<OneD, NekDouble> &outarray)
197{
198 switch (dir)
199 {
200 case 0:
201 {
202 v_PhysDeriv(inarray, outarray, NullNekDouble1DArray,
204 break;
205 }
206 case 1:
207 {
208 v_PhysDeriv(inarray, NullNekDouble1DArray, outarray,
210 break;
211 }
212 case 2:
213 {
215 outarray);
216 break;
217 }
218 default:
219 {
220 ASSERTL1(false, "input dir is out of range");
221 }
222 break;
223 }
224}
225
230{
231 StdTetExp::v_PhysDeriv(inarray, out_d0, out_d1, out_d2);
232}
233
234//---------------------------------------
235// Transforms
236//---------------------------------------
237
238/**
239 * @note 'r' (base[2]) runs fastest in this element
240 *
241 * \f$ u^{\delta} (\xi_{1i}, \xi_{2j}, \xi_{3k}) = \sum_{m(pqr)} \hat
242 * u_{pqr} \phi_{pqr} (\xi_{1i}, \xi_{2j}, \xi_{3k})\f$
243 *
244 * Backward transformation is three dimensional tensorial expansion
245 * \f$ u (\xi_{1i}, \xi_{2j}, \xi_{3k}) = \sum_{p=0}^{Q_x} \psi_p^a
246 * (\xi_{1i}) \lbrace { \sum_{q=0}^{Q_y} \psi_{pq}^b (\xi_{2j})
247 * \lbrace { \sum_{r=0}^{Q_z} \hat u_{pqr} \psi_{pqr}^c (\xi_{3k})
248 * \rbrace} \rbrace}. \f$ And sumfactorizing step of the form is as:\\
249 *
250 * \f$ f_{pq} (\xi_{3k}) = \sum_{r=0}^{Q_z} \hat u_{pqr} \psi_{pqr}^c
251 * (\xi_{3k}),\\
252 *
253 * g_{p} (\xi_{2j}, \xi_{3k}) = \sum_{r=0}^{Q_y} \psi_{pq}^b
254 * (\xi_{2j}) f_{pq} (\xi_{3k}),\\
255 *
256 * u(\xi_{1i}, \xi_{2j}, \xi_{3k}) = \sum_{p=0}^{Q_x} \psi_{p}^a
257 * (\xi_{1i}) g_{p} (\xi_{2j}, \xi_{3k}). \f$
258 */
260 Array<OneD, NekDouble> &outarray)
261{
264 "Basis[1] is not a general tensor type");
265
268 "Basis[2] is not a general tensor type");
269
270 if (m_base[0]->Collocation() && m_base[1]->Collocation() &&
271 m_base[2]->Collocation())
272 {
274 m_base[2]->GetNumPoints(),
275 inarray, 1, outarray, 1);
276 }
277 else
278 {
279 StdTetExp::v_BwdTrans_SumFac(inarray, outarray);
280 }
281}
282
283/**
284 * Sum-factorisation implementation of the BwdTrans operation.
285 */
287 Array<OneD, NekDouble> &outarray)
288{
289 int nquad1 = m_base[1]->GetNumPoints();
290 int nquad2 = m_base[2]->GetNumPoints();
291 int order0 = m_base[0]->GetNumModes();
292 int order1 = m_base[1]->GetNumModes();
293
294 Array<OneD, NekDouble> wsp(nquad2 * order0 * (2 * order1 - order0 + 1) / 2 +
295 nquad2 * nquad1 * order0);
296
297 BwdTrans_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetBdata(),
298 m_base[2]->GetBdata(), inarray, outarray, wsp, true,
299 true, true);
300}
301
302/**
303 * @param base0 x-dirn basis matrix
304 * @param base1 y-dirn basis matrix
305 * @param base2 z-dirn basis matrix
306 * @param inarray Input vector of modes.
307 * @param outarray Output vector of physical space data.
308 * @param wsp Workspace of size Q_x*P_z*(P_y+Q_y)
309 * @param doCheckCollDir0 Check for collocation of basis.
310 * @param doCheckCollDir1 Check for collocation of basis.
311 * @param doCheckCollDir2 Check for collocation of basis.
312 * @todo Account for some directions being collocated. See
313 * StdQuadExp as an example.
314 */
316 const Array<OneD, const NekDouble> &base0,
317 const Array<OneD, const NekDouble> &base1,
318 const Array<OneD, const NekDouble> &base2,
319 const Array<OneD, const NekDouble> &inarray,
321 [[maybe_unused]] bool doCheckCollDir0,
322 [[maybe_unused]] bool doCheckCollDir1,
323 [[maybe_unused]] bool doCheckCollDir2)
324{
325 int nquad0 = m_base[0]->GetNumPoints();
326 int nquad1 = m_base[1]->GetNumPoints();
327 int nquad2 = m_base[2]->GetNumPoints();
328
329 int order0 = m_base[0]->GetNumModes();
330 int order1 = m_base[1]->GetNumModes();
331 int order2 = m_base[2]->GetNumModes();
332
333 Array<OneD, NekDouble> tmp = wsp;
335 tmp + nquad2 * order0 * (2 * order1 - order0 + 1) / 2;
336
337 int i, j, mode, mode1, cnt;
338
339 // Perform summation over '2' direction
340 mode = mode1 = cnt = 0;
341 for (i = 0; i < order0; ++i)
342 {
343 for (j = 0; j < order1 - i; ++j, ++cnt)
344 {
345 Blas::Dgemv('N', nquad2, order2 - i - j, 1.0,
346 base2.data() + mode * nquad2, nquad2,
347 inarray.data() + mode1, 1, 0.0,
348 tmp.data() + cnt * nquad2, 1);
349 mode += order2 - i - j;
350 mode1 += order2 - i - j;
351 }
352 // increment mode in case order1!=order2
353 for (j = order1 - i; j < order2 - i; ++j)
354 {
355 mode += order2 - i - j;
356 }
357 }
358
359 // fix for modified basis by adding split of top singular
360 // vertex mode - currently (1+c)/2 x (1-b)/2 x (1-a)/2
361 // component is evaluated
363 {
364 // top singular vertex - (1+c)/2 x (1+b)/2 x (1-a)/2 component
365 Blas::Daxpy(nquad2, inarray[1], base2.data() + nquad2, 1,
366 &tmp[0] + nquad2, 1);
367
368 // top singular vertex - (1+c)/2 x (1-b)/2 x (1+a)/2 component
369 Blas::Daxpy(nquad2, inarray[1], base2.data() + nquad2, 1,
370 &tmp[0] + order1 * nquad2, 1);
371 }
372
373 // Perform summation over '1' direction
374 mode = 0;
375 for (i = 0; i < order0; ++i)
376 {
377 Blas::Dgemm('N', 'T', nquad1, nquad2, order1 - i, 1.0,
378 base1.data() + mode * nquad1, nquad1,
379 tmp.data() + mode * nquad2, nquad2, 0.0,
380 tmp1.data() + i * nquad1 * nquad2, nquad1);
381 mode += order1 - i;
382 }
383
384 // fix for modified basis by adding additional split of
385 // top and base singular vertex modes as well as singular
386 // edge
388 {
389 // use tmp to sort out singular vertices and
390 // singular edge components with (1+b)/2 (1+a)/2 form
391 for (i = 0; i < nquad2; ++i)
392 {
393 Blas::Daxpy(nquad1, tmp[nquad2 + i], base1.data() + nquad1, 1,
394 &tmp1[nquad1 * nquad2] + i * nquad1, 1);
395 }
396 }
397
398 // Perform summation over '0' direction
399 Blas::Dgemm('N', 'T', nquad0, nquad1 * nquad2, order0, 1.0, base0.data(),
400 nquad0, tmp1.data(), nquad1 * nquad2, 0.0, outarray.data(),
401 nquad0);
402}
403
404/**
405 * @param inarray array of physical quadrature points to be
406 * transformed.
407 * @param outarray updated array of expansion coefficients.
408 */
410 Array<OneD, NekDouble> &outarray)
411{ // int numMax = nmodes0;
412 v_IProductWRTBase(inarray, outarray);
413
414 // get Mass matrix inverse
415 StdMatrixKey masskey(eInvMass, DetShapeType(), *this);
416 DNekMatSharedPtr matsys = GetStdMatrix(masskey);
417
418 // copy inarray in case inarray == outarray
419 DNekVec in(m_ncoeffs, outarray);
420 DNekVec out(m_ncoeffs, outarray, eWrapper);
421
422 out = (*matsys) * in;
423}
424
425//---------------------------------------
426// Inner product functions
427//---------------------------------------
428
429/**
430 * \f$ \begin{array}{rcl} I_{pqr} = (\phi_{pqr}, u)_{\delta} & = &
431 * \sum_{i=0}^{nq_0} \sum_{j=0}^{nq_1} \sum_{k=0}^{nq_2} \psi_{p}^{a}
432 * (\eta_{1i}) \psi_{pq}^{b} (\eta_{2j}) \psi_{pqr}^{c} (\eta_{3k})
433 * w_i w_j w_k u(\eta_{1,i} \eta_{2,j} \eta_{3,k}) J_{i,j,k}\\ & = &
434 * \sum_{i=0}^{nq_0} \psi_p^a(\eta_{1,i}) \sum_{j=0}^{nq_1}
435 * \psi_{pq}^b(\eta_{2,j}) \sum_{k=0}^{nq_2} \psi_{pqr}^c
436 * u(\eta_{1i},\eta_{2j},\eta_{3k}) J_{i,j,k} \end{array} \f$ \n
437 *
438 * where
439 *
440 * \f$ \phi_{pqr} (\xi_1 , \xi_2 , \xi_3) = \psi_p^a (\eta_1)
441 * \psi_{pq}^b (\eta_2) \psi_{pqr}^c (\eta_3) \f$
442 *
443 * which can be implemented as \n \f$f_{pqr} (\xi_{3k}) =
444 * \sum_{k=0}^{nq_3} \psi_{pqr}^c u(\eta_{1i},\eta_{2j},\eta_{3k})
445 *
446 * J_{i,j,k} = {\bf B_3 U} \f$ \n
447 *
448 * \f$ g_{pq} (\xi_{3k}) = \sum_{j=0}^{nq_1} \psi_{pq}^b (\xi_{2j})
449 * f_{pqr} (\xi_{3k}) = {\bf B_2 F} \f$ \n
450 *
451 * \f$ (\phi_{pqr}, u)_{\delta} = \sum_{k=0}^{nq_0} \psi_{p}^a
452 * (\xi_{3k}) g_{pq} (\xi_{3k}) = {\bf B_1 G} \f$
453 *
454 * @param inarray Function evaluated at physical collocation
455 * points.
456 * @param outarray Inner product with respect to each basis
457 * function over the element.
458 */
460 Array<OneD, NekDouble> &outarray)
461{
464 "Basis[1] is not a general tensor type");
465
468 "Basis[2] is not a general tensor type");
469
470 if (m_base[0]->Collocation() && m_base[1]->Collocation())
471 {
472 MultiplyByQuadratureMetric(inarray, outarray);
473 }
474 else
475 {
476 StdTetExp::v_IProductWRTBase_SumFac(inarray, outarray);
477 }
478}
479
480/**
481 * @param inarray Function evaluated at physical collocation
482 * points.
483 * @param outarray Inner product with respect to each basis
484 * function over the element.
485 */
487 const Array<OneD, const NekDouble> &inarray,
488 Array<OneD, NekDouble> &outarray, bool multiplybyweights)
489{
490 int nquad0 = m_base[0]->GetNumPoints();
491 int nquad1 = m_base[1]->GetNumPoints();
492 int nquad2 = m_base[2]->GetNumPoints();
493 int order0 = m_base[0]->GetNumModes();
494 int order1 = m_base[1]->GetNumModes();
495
496 Array<OneD, NekDouble> wsp(nquad1 * nquad2 * order0 +
497 nquad2 * order0 * (2 * order1 - order0 + 1) / 2);
498
499 if (multiplybyweights)
500 {
501 Array<OneD, NekDouble> tmp(nquad0 * nquad1 * nquad2);
502 MultiplyByQuadratureMetric(inarray, tmp);
503
505 m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
506 tmp, outarray, wsp, true, true, true);
507 }
508 else
509 {
511 m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
512 inarray, outarray, wsp, true, true, true);
513 }
514}
515
517 const Array<OneD, const NekDouble> &base0,
518 const Array<OneD, const NekDouble> &base1,
519 const Array<OneD, const NekDouble> &base2,
520 const Array<OneD, const NekDouble> &inarray,
522 [[maybe_unused]] bool doCheckCollDir0,
523 [[maybe_unused]] bool doCheckCollDir1,
524 [[maybe_unused]] bool doCheckCollDir2)
525{
526 int nquad0 = m_base[0]->GetNumPoints();
527 int nquad1 = m_base[1]->GetNumPoints();
528 int nquad2 = m_base[2]->GetNumPoints();
529
530 int order0 = m_base[0]->GetNumModes();
531 int order1 = m_base[1]->GetNumModes();
532 int order2 = m_base[2]->GetNumModes();
533
534 Array<OneD, NekDouble> tmp1 = wsp;
535 Array<OneD, NekDouble> tmp2 = wsp + nquad1 * nquad2 * order0;
536
537 int i, j, mode, mode1, cnt;
538
539 // Inner product with respect to the '0' direction
540 Blas::Dgemm('T', 'N', nquad1 * nquad2, order0, nquad0, 1.0, inarray.data(),
541 nquad0, base0.data(), nquad0, 0.0, tmp1.data(),
542 nquad1 * nquad2);
543
544 // Inner product with respect to the '1' direction
545 for (mode = i = 0; i < order0; ++i)
546 {
547 Blas::Dgemm('T', 'N', nquad2, order1 - i, nquad1, 1.0,
548 tmp1.data() + i * nquad1 * nquad2, nquad1,
549 base1.data() + mode * nquad1, nquad1, 0.0,
550 tmp2.data() + mode * nquad2, nquad2);
551 mode += order1 - i;
552 }
553
554 // fix for modified basis for base singular vertex
556 {
557 // base singular vertex and singular edge (1+b)/2
558 //(1+a)/2 components (makes tmp[nquad2] entry into (1+b)/2)
559 Blas::Dgemv('T', nquad1, nquad2, 1.0, tmp1.data() + nquad1 * nquad2,
560 nquad1, base1.data() + nquad1, 1, 1.0, tmp2.data() + nquad2,
561 1);
562 }
563
564 // Inner product with respect to the '2' direction
565 mode = mode1 = cnt = 0;
566 for (i = 0; i < order0; ++i)
567 {
568 for (j = 0; j < order1 - i; ++j, ++cnt)
569 {
570 Blas::Dgemv('T', nquad2, order2 - i - j, 1.0,
571 base2.data() + mode * nquad2, nquad2,
572 tmp2.data() + cnt * nquad2, 1, 0.0,
573 outarray.data() + mode1, 1);
574 mode += order2 - i - j;
575 mode1 += order2 - i - j;
576 }
577 // increment mode in case order1!=order2
578 for (j = order1 - i; j < order2 - i; ++j)
579 {
580 mode += order2 - i - j;
581 }
582 }
583
584 // fix for modified basis for top singular vertex component
585 // Already have evaluated (1+c)/2 (1-b)/2 (1-a)/2
587 {
588 // add in (1+c)/2 (1+b)/2 component
589 outarray[1] +=
590 Blas::Ddot(nquad2, base2.data() + nquad2, 1, &tmp2[nquad2], 1);
591
592 // add in (1+c)/2 (1-b)/2 (1+a)/2 component
593 outarray[1] += Blas::Ddot(nquad2, base2.data() + nquad2, 1,
594 &tmp2[nquad2 * order1], 1);
595 }
596}
597
599 const int dir, const Array<OneD, const NekDouble> &inarray,
600 Array<OneD, NekDouble> &outarray)
601{
602 StdTetExp::v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
603}
604
605/**
606 * @param inarray Function evaluated at physical collocation
607 * points.
608 * @param outarray Inner product with respect to each basis
609 * function over the element.
610 */
612 const int dir, const Array<OneD, const NekDouble> &inarray,
613 Array<OneD, NekDouble> &outarray)
614{
615 int i;
616 int nquad0 = m_base[0]->GetNumPoints();
617 int nquad1 = m_base[1]->GetNumPoints();
618 int nquad2 = m_base[2]->GetNumPoints();
619 int nqtot = nquad0 * nquad1 * nquad2;
620 int nmodes0 = m_base[0]->GetNumModes();
621 int nmodes1 = m_base[1]->GetNumModes();
622 int wspsize = nquad0 + nquad1 + nquad2 + max(nqtot, m_ncoeffs) +
623 nquad1 * nquad2 * nmodes0 +
624 nquad2 * nmodes0 * (2 * nmodes1 - nmodes0 + 1) / 2;
625
626 Array<OneD, NekDouble> gfac0(wspsize);
627 Array<OneD, NekDouble> gfac1(gfac0 + nquad0);
628 Array<OneD, NekDouble> gfac2(gfac1 + nquad1);
629 Array<OneD, NekDouble> tmp0(gfac2 + nquad2);
630 Array<OneD, NekDouble> wsp(tmp0 + max(nqtot, m_ncoeffs));
631
632 const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
633 const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
634 const Array<OneD, const NekDouble> &z2 = m_base[2]->GetZ();
635
636 // set up geometric factor: (1+z0)/2
637 for (i = 0; i < nquad0; ++i)
638 {
639 gfac0[i] = 0.5 * (1 + z0[i]);
640 }
641
642 // set up geometric factor: 2/(1-z1)
643 for (i = 0; i < nquad1; ++i)
644 {
645 gfac1[i] = 2.0 / (1 - z1[i]);
646 }
647
648 // Set up geometric factor: 2/(1-z2)
649 for (i = 0; i < nquad2; ++i)
650 {
651 gfac2[i] = 2.0 / (1 - z2[i]);
652 }
653
654 // Derivative in first direction is always scaled as follows
655 for (i = 0; i < nquad1 * nquad2; ++i)
656 {
657 Vmath::Smul(nquad0, gfac1[i % nquad1], &inarray[0] + i * nquad0, 1,
658 &tmp0[0] + i * nquad0, 1);
659 }
660 for (i = 0; i < nquad2; ++i)
661 {
662 Vmath::Smul(nquad0 * nquad1, gfac2[i], &tmp0[0] + i * nquad0 * nquad1,
663 1, &tmp0[0] + i * nquad0 * nquad1, 1);
664 }
665
666 MultiplyByQuadratureMetric(tmp0, tmp0);
667
668 switch (dir)
669 {
670 case 0:
671 {
673 m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
674 m_base[2]->GetBdata(), tmp0, outarray, wsp, false, true, true);
675 }
676 break;
677 case 1:
678 {
680
681 for (i = 0; i < nquad1 * nquad2; ++i)
682 {
683 Vmath::Vmul(nquad0, &gfac0[0], 1, &tmp0[0] + i * nquad0, 1,
684 &tmp0[0] + i * nquad0, 1);
685 }
686
688 m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
689 m_base[2]->GetBdata(), tmp0, tmp3, wsp, false, true, true);
690
691 for (i = 0; i < nquad2; ++i)
692 {
693 Vmath::Smul(nquad0 * nquad1, gfac2[i],
694 &inarray[0] + i * nquad0 * nquad1, 1,
695 &tmp0[0] + i * nquad0 * nquad1, 1);
696 }
697 MultiplyByQuadratureMetric(tmp0, tmp0);
699 m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
700 m_base[2]->GetBdata(), tmp0, outarray, wsp, true, false, true);
701 Vmath::Vadd(m_ncoeffs, &tmp3[0], 1, &outarray[0], 1, &outarray[0],
702 1);
703 }
704 break;
705 case 2:
706 {
709 for (i = 0; i < nquad1; ++i)
710 {
711 gfac1[i] = (1 + z1[i]) / 2;
712 }
713
714 for (i = 0; i < nquad1 * nquad2; ++i)
715 {
716 Vmath::Vmul(nquad0, &gfac0[0], 1, &tmp0[0] + i * nquad0, 1,
717 &tmp0[0] + i * nquad0, 1);
718 }
720 m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
721 m_base[2]->GetBdata(), tmp0, tmp3, wsp, false, true, true);
722
723 for (i = 0; i < nquad2; ++i)
724 {
725 Vmath::Smul(nquad0 * nquad1, gfac2[i],
726 &inarray[0] + i * nquad0 * nquad1, 1,
727 &tmp0[0] + i * nquad0 * nquad1, 1);
728 }
729 for (i = 0; i < nquad1 * nquad2; ++i)
730 {
731 Vmath::Smul(nquad0, gfac1[i % nquad1], &tmp0[0] + i * nquad0, 1,
732 &tmp0[0] + i * nquad0, 1);
733 }
734 MultiplyByQuadratureMetric(tmp0, tmp0);
736 m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
737 m_base[2]->GetBdata(), tmp0, tmp4, wsp, true, false, true);
738
739 MultiplyByQuadratureMetric(inarray, tmp0);
741 m_base[0]->GetBdata(), m_base[1]->GetBdata(),
742 m_base[2]->GetDbdata(), tmp0, outarray, wsp, true, true, false);
743
744 Vmath::Vadd(m_ncoeffs, &tmp3[0], 1, &outarray[0], 1, &outarray[0],
745 1);
746 Vmath::Vadd(m_ncoeffs, &tmp4[0], 1, &outarray[0], 1, &outarray[0],
747 1);
748 }
749 break;
750 default:
751 {
752 ASSERTL1(false, "input dir is out of range");
753 }
754 break;
755 }
756}
757
758//---------------------------------------
759// Evaluation functions
760//---------------------------------------
761
764{
765 NekDouble d2 = 1.0 - xi[2];
766 NekDouble d12 = -xi[1] - xi[2];
767 if (fabs(d2) < NekConstants::kNekZeroTol)
768 {
769 if (d2 >= 0.)
770 {
772 }
773 else
774 {
776 }
777 }
778 if (fabs(d12) < NekConstants::kNekZeroTol)
779 {
780 if (d12 >= 0.)
781 {
783 }
784 else
785 {
787 }
788 }
789 eta[0] = 2.0 * (1.0 + xi[0]) / d12 - 1.0;
790 eta[1] = 2.0 * (1.0 + xi[1]) / d2 - 1.0;
791 eta[2] = xi[2];
792}
793
796{
797 xi[2] = eta[2];
798 xi[1] = (1.0 + eta[1]) * (1.0 - xi[2]) * 0.5 - 1.0;
799 xi[0] = (1.0 + eta[0]) * (-xi[1] - xi[2]) * 0.5 - 1.0;
800}
801
802void StdTetExp::v_FillMode(const int mode, Array<OneD, NekDouble> &outarray)
803{
805 tmp[mode] = 1.0;
806 StdTetExp::v_BwdTrans(tmp, outarray);
807}
808
810 const Array<OneD, const NekDouble> &coords, int mode)
811{
813 LocCoordToLocCollapsed(coords, coll);
814
815 const int nm1 = m_base[1]->GetNumModes();
816 const int nm2 = m_base[2]->GetNumModes();
817
818 const int b = 2 * nm2 + 1;
819 const int mode0 = floor(0.5 * (b - sqrt(b * b - 8.0 * mode / nm1)));
820 const int tmp =
821 mode - nm1 * (mode0 * (nm2 - 1) + 1 - (mode0 - 2) * (mode0 - 1) / 2);
822 const int mode1 = tmp / (nm2 - mode0);
823 const int mode2 = tmp % (nm2 - mode0);
824
826 {
827 // Handle the collapsed vertices and edges in the modified
828 // basis.
829 if (mode == 1)
830 {
831 // Collapsed top vertex
832 return StdExpansion::BaryEvaluateBasis<2>(coll[2], 1);
833 }
834 else if (mode0 == 0 && mode2 == 1)
835 {
836 return StdExpansion::BaryEvaluateBasis<1>(coll[1], 0) *
837 StdExpansion::BaryEvaluateBasis<2>(coll[2], 1);
838 }
839 else if (mode0 == 1 && mode1 == 1 && mode2 == 0)
840 {
841 return StdExpansion::BaryEvaluateBasis<0>(coll[0], 0) *
842 StdExpansion::BaryEvaluateBasis<1>(coll[1], 1);
843 }
844 }
845
846 return StdExpansion::BaryEvaluateBasis<0>(coll[0], mode0) *
847 StdExpansion::BaryEvaluateBasis<1>(coll[1], mode1) *
848 StdExpansion::BaryEvaluateBasis<2>(coll[2], mode2);
849}
850
852 const Array<OneD, NekDouble> &coord,
853 const Array<OneD, const NekDouble> &inarray,
854 std::array<NekDouble, 3> &firstOrderDerivs)
855{
856 // Collapse coordinates
857 Array<OneD, NekDouble> coll(3, 0.0);
858 LocCoordToLocCollapsed(coord, coll);
859
860 // If near singularity do the old interpolation matrix method
861 if ((1 - coll[1]) < 1e-5 || (1 - coll[2]) < 1e-5)
862 {
863 int totPoints = GetTotPoints();
864 Array<OneD, NekDouble> EphysDeriv0(totPoints), EphysDeriv1(totPoints),
865 EphysDeriv2(totPoints);
866 PhysDeriv(inarray, EphysDeriv0, EphysDeriv1, EphysDeriv2);
867
869 I[0] = GetBase()[0]->GetI(coll);
870 I[1] = GetBase()[1]->GetI(coll + 1);
871 I[2] = GetBase()[2]->GetI(coll + 2);
872
873 firstOrderDerivs[0] = PhysEvaluate(I, EphysDeriv0);
874 firstOrderDerivs[1] = PhysEvaluate(I, EphysDeriv1);
875 firstOrderDerivs[2] = PhysEvaluate(I, EphysDeriv2);
876 return PhysEvaluate(I, inarray);
877 }
878
879 std::array<NekDouble, 3> interDeriv;
880 NekDouble val = BaryTensorDeriv(coll, inarray, interDeriv);
881
882 // calculate 2.0/((1-eta_1)(1-eta_2)) * Out_dEta0
883 NekDouble temp = 2.0 / ((1 - coll[1]) * (1 - coll[2]));
884 interDeriv[0] *= temp;
885
886 // out_dxi0 = 4.0/((1-eta_1)(1-eta_2)) * Out_dEta0
887 firstOrderDerivs[0] = 2 * interDeriv[0];
888
889 // fac0 = 1 + eta_0
890 NekDouble fac0;
891 fac0 = 1 + coll[0];
892
893 // calculate 2.0*(1+eta_0)/((1-eta_1)(1-eta_2)) * Out_dEta0
894 interDeriv[0] *= fac0;
895
896 // calculate 2/(1.0-eta_2) * out_dEta1
897 fac0 = 2 / (1 - coll[2]);
898 interDeriv[1] *= fac0;
899
900 // calculate out_dxi1 = 2.0(1+eta_0)/((1-eta_1)(1-eta_2))
901 // * Out_dEta0 + 2/(1.0-eta_2) out_dEta1
902 firstOrderDerivs[1] = interDeriv[0] + interDeriv[1];
903
904 // calculate (1 + eta_1)/(1 -eta_2)*out_dEta1
905 fac0 = (1 + coll[1]) / 2;
906 interDeriv[1] *= fac0;
907
908 // calculate out_dxi2 =
909 // 2.0(1+eta_0)/((1-eta_1)(1-eta_2)) Out_dEta0 +
910 // (1 + eta_1)/(1 -eta_2)*out_dEta1 + out_dEta2
911 firstOrderDerivs[2] = interDeriv[0] + interDeriv[1] + interDeriv[2];
912
913 return val;
914}
915
916void StdTetExp::v_GetTraceNumModes(const int fid, int &numModes0,
917 int &numModes1,
918 [[maybe_unused]] Orientation faceOrient)
919{
920 int nummodes[3] = {m_base[0]->GetNumModes(), m_base[1]->GetNumModes(),
921 m_base[2]->GetNumModes()};
922 switch (fid)
923 {
924 case 0:
925 {
926 numModes0 = nummodes[0];
927 numModes1 = nummodes[1];
928 }
929 break;
930 case 1:
931 {
932 numModes0 = nummodes[0];
933 numModes1 = nummodes[2];
934 }
935 break;
936 case 2:
937 case 3:
938 {
939 numModes0 = nummodes[1];
940 numModes1 = nummodes[2];
941 }
942 break;
943 }
944}
945
946//---------------------------
947// Helper functions
948//---------------------------
949
951{
952 return 4;
953}
954
956{
957 return 6;
958}
959
961{
962 return 4;
963}
964
969
971{
974 "BasisType is not a boundary interior form");
977 "BasisType is not a boundary interior form");
980 "BasisType is not a boundary interior form");
981
982 int P = m_base[0]->GetNumModes();
983 int Q = m_base[1]->GetNumModes();
984 int R = m_base[2]->GetNumModes();
985
987}
988
990{
993 "BasisType is not a boundary interior form");
996 "BasisType is not a boundary interior form");
999 "BasisType is not a boundary interior form");
1000
1001 int P = m_base[0]->GetNumModes() - 1;
1002 int Q = m_base[1]->GetNumModes() - 1;
1003 int R = m_base[2]->GetNumModes() - 1;
1004
1005 return (Q + 1) + P * (1 + 2 * Q - P) / 2 // base face
1006 + (R + 1) + P * (1 + 2 * R - P) / 2 // front face
1007 + 2 * (R + 1) + Q * (1 + 2 * R - Q); // back two faces
1008}
1009
1010int StdTetExp::v_GetTraceNcoeffs(const int i) const
1011{
1012 ASSERTL2((i >= 0) && (i <= 3), "face id is out of range");
1013 int nFaceCoeffs = 0;
1014 int nummodesA, nummodesB, P, Q;
1015 if (i == 0)
1016 {
1017 nummodesA = GetBasisNumModes(0);
1018 nummodesB = GetBasisNumModes(1);
1019 }
1020 else if ((i == 1) || (i == 2))
1021 {
1022 nummodesA = GetBasisNumModes(0);
1023 nummodesB = GetBasisNumModes(2);
1024 }
1025 else
1026 {
1027 nummodesA = GetBasisNumModes(1);
1028 nummodesB = GetBasisNumModes(2);
1029 }
1030 P = nummodesA - 1;
1031 Q = nummodesB - 1;
1032 nFaceCoeffs = Q + 1 + (P * (1 + 2 * Q - P)) / 2;
1033 return nFaceCoeffs;
1034}
1035
1037{
1038 ASSERTL2((i >= 0) && (i <= 3), "face id is out of range");
1039 int Pi = m_base[0]->GetNumModes() - 2;
1040 int Qi = m_base[1]->GetNumModes() - 2;
1041 int Ri = m_base[2]->GetNumModes() - 2;
1042
1043 if ((i == 0))
1044 {
1045 return Pi * (2 * Qi - Pi - 1) / 2;
1046 }
1047 else if ((i == 1))
1048 {
1049 return Pi * (2 * Ri - Pi - 1) / 2;
1050 }
1051 else
1052 {
1053 return Qi * (2 * Ri - Qi - 1) / 2;
1054 }
1055}
1056
1058{
1059 ASSERTL2(i >= 0 && i <= 3, "face id is out of range");
1060
1061 if (i == 0)
1062 {
1063 return m_base[0]->GetNumPoints() * m_base[1]->GetNumPoints();
1064 }
1065 else if (i == 1)
1066 {
1067 return m_base[0]->GetNumPoints() * m_base[2]->GetNumPoints();
1068 }
1069 else
1070 {
1071 return m_base[1]->GetNumPoints() * m_base[2]->GetNumPoints();
1072 }
1073}
1074
1075int StdTetExp::v_GetEdgeNcoeffs(const int i) const
1076{
1077 ASSERTL2((i >= 0) && (i <= 5), "edge id is out of range");
1078 int P = m_base[0]->GetNumModes();
1079 int Q = m_base[1]->GetNumModes();
1080 int R = m_base[2]->GetNumModes();
1081
1082 if (i == 0)
1083 {
1084 return P;
1085 }
1086 else if (i == 1 || i == 2)
1087 {
1088 return Q;
1089 }
1090 else
1091 {
1092 return R;
1093 }
1094}
1095
1097 const int j) const
1098{
1099 ASSERTL2(i >= 0 && i <= 3, "face id is out of range");
1100 ASSERTL2(j == 0 || j == 1, "face direction is out of range");
1101
1102 if (i == 0)
1103 {
1104 return m_base[j]->GetPointsKey();
1105 }
1106 else if (i == 1)
1107 {
1108 return m_base[2 * j]->GetPointsKey();
1109 }
1110 else
1111 {
1112 return m_base[j + 1]->GetPointsKey();
1113 }
1114}
1115
1117 const std::vector<unsigned int> &nummodes, int &modes_offset)
1118{
1120 nummodes[modes_offset], nummodes[modes_offset + 1],
1121 nummodes[modes_offset + 2]);
1122 modes_offset += 3;
1123
1124 return nmodes;
1125}
1126
1128 const int k,
1129 bool UseGLL) const
1130{
1131 ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
1132 ASSERTL2(k == 0 || k == 1, "face direction out of range");
1133
1134 int dir = k;
1135 switch (i)
1136 {
1137 case 0:
1138 dir = k;
1139 break;
1140 case 1:
1141 dir = 2 * k;
1142 break;
1143 case 2:
1144 case 3:
1145 dir = k + 1;
1146 break;
1147 }
1148
1149 return EvaluateTriFaceBasisKey(k, m_base[dir], UseGLL);
1150}
1151
1155{
1156 Array<OneD, const NekDouble> eta_x = m_base[0]->GetZ();
1157 Array<OneD, const NekDouble> eta_y = m_base[1]->GetZ();
1158 Array<OneD, const NekDouble> eta_z = m_base[2]->GetZ();
1159 int Qx = GetNumPoints(0);
1160 int Qy = GetNumPoints(1);
1161 int Qz = GetNumPoints(2);
1162
1163 // Convert collapsed coordinates into cartesian coordinates: eta
1164 // --> xi
1165 for (int k = 0; k < Qz; ++k)
1166 {
1167 for (int j = 0; j < Qy; ++j)
1168 {
1169 for (int i = 0; i < Qx; ++i)
1170 {
1171 int s = i + Qx * (j + Qy * k);
1172 xi_x[s] =
1173 (eta_x[i] + 1.0) * (1.0 - eta_y[j]) * (1.0 - eta_z[k]) / 4 -
1174 1.0;
1175 xi_y[s] = (eta_y[j] + 1.0) * (1.0 - eta_z[k]) / 2 - 1.0;
1176 xi_z[s] = eta_z[k];
1177 }
1178 }
1179 }
1180}
1181
1183{
1184 return (m_base[0]->GetBasisType() == LibUtilities::eModified_A) &&
1185 (m_base[1]->GetBasisType() == LibUtilities::eModified_B) &&
1187}
1188
1189//--------------------------
1190// Mappings
1191//--------------------------
1192int StdTetExp::v_GetVertexMap(const int localVertexId, bool useCoeffPacking)
1193{
1197 "Mapping not defined for this type of basis");
1198
1199 int localDOF = 0;
1200 if (useCoeffPacking == true) // follow packing of coefficients i.e q,r,p
1201 {
1202 switch (localVertexId)
1203 {
1204 case 0:
1205 {
1206 localDOF = GetMode(0, 0, 0);
1207 break;
1208 }
1209 case 1:
1210 {
1211 localDOF = GetMode(0, 0, 1);
1212 break;
1213 }
1214 case 2:
1215 {
1216 localDOF = GetMode(0, 1, 0);
1217 break;
1218 }
1219 case 3:
1220 {
1221 localDOF = GetMode(1, 0, 0);
1222 break;
1223 }
1224 default:
1225 {
1226 ASSERTL0(false, "Vertex ID must be between 0 and 3");
1227 break;
1228 }
1229 }
1230 }
1231 else
1232 {
1233 switch (localVertexId)
1234 {
1235 case 0:
1236 {
1237 localDOF = GetMode(0, 0, 0);
1238 break;
1239 }
1240 case 1:
1241 {
1242 localDOF = GetMode(1, 0, 0);
1243 break;
1244 }
1245 case 2:
1246 {
1247 localDOF = GetMode(0, 1, 0);
1248 break;
1249 }
1250 case 3:
1251 {
1252 localDOF = GetMode(0, 0, 1);
1253 break;
1254 }
1255 default:
1256 {
1257 ASSERTL0(false, "Vertex ID must be between 0 and 3");
1258 break;
1259 }
1260 }
1261 }
1262
1263 return localDOF;
1264}
1265
1266/**
1267 * Maps interior modes of an edge to the elemental modes.
1268 */
1269
1270/**
1271 * List of all interior modes in the expansion.
1272 */
1274{
1277 "BasisType is not a boundary interior form");
1280 "BasisType is not a boundary interior form");
1283 "BasisType is not a boundary interior form");
1284
1285 int P = m_base[0]->GetNumModes();
1286 int Q = m_base[1]->GetNumModes();
1287 int R = m_base[2]->GetNumModes();
1288
1289 int nIntCoeffs = m_ncoeffs - NumBndryCoeffs();
1290
1291 if (outarray.size() != nIntCoeffs)
1292 {
1293 outarray = Array<OneD, unsigned int>(nIntCoeffs);
1294 }
1295
1296 int idx = 0;
1297 for (int i = 2; i < P; ++i)
1298 {
1299 for (int j = 1; j < Q - i - 1; ++j)
1300 {
1301 for (int k = 1; k < R - i - j; ++k)
1302 {
1303 outarray[idx++] = GetMode(i, j, k);
1304 }
1305 }
1306 }
1307}
1308
1309/**
1310 * List of all boundary modes in the the expansion.
1311 */
1313{
1316 "BasisType is not a boundary interior form");
1319 "BasisType is not a boundary interior form");
1322 "BasisType is not a boundary interior form");
1323
1324 int P = m_base[0]->GetNumModes();
1325 int Q = m_base[1]->GetNumModes();
1326 int R = m_base[2]->GetNumModes();
1327
1328 int i, j, k;
1329 int idx = 0;
1330
1331 int nBnd = NumBndryCoeffs();
1332
1333 if (outarray.size() != nBnd)
1334 {
1335 outarray = Array<OneD, unsigned int>(nBnd);
1336 }
1337
1338 for (i = 0; i < P; ++i)
1339 {
1340 // First two Q-R planes are entirely boundary modes
1341 if (i < 2)
1342 {
1343 for (j = 0; j < Q - i; j++)
1344 {
1345 for (k = 0; k < R - i - j; ++k)
1346 {
1347 outarray[idx++] = GetMode(i, j, k);
1348 }
1349 }
1350 }
1351 // Remaining Q-R planes contain boundary modes on bottom and
1352 // left edge.
1353 else
1354 {
1355 for (k = 0; k < R - i; ++k)
1356 {
1357 outarray[idx++] = GetMode(i, 0, k);
1358 }
1359 for (j = 1; j < Q - i; ++j)
1360 {
1361 outarray[idx++] = GetMode(i, j, 0);
1362 }
1363 }
1364 }
1365}
1366
1367void StdTetExp::v_GetTraceCoeffMap(const unsigned int fid,
1368 Array<OneD, unsigned int> &maparray)
1369{
1370 int i, j, k;
1371 int P = 0, Q = 0, idx = 0;
1372 int nFaceCoeffs = 0;
1373
1374 switch (fid)
1375 {
1376 case 0:
1377 P = m_base[0]->GetNumModes();
1378 Q = m_base[1]->GetNumModes();
1379 break;
1380 case 1:
1381 P = m_base[0]->GetNumModes();
1382 Q = m_base[2]->GetNumModes();
1383 break;
1384 case 2:
1385 case 3:
1386 P = m_base[1]->GetNumModes();
1387 Q = m_base[2]->GetNumModes();
1388 break;
1389 default:
1390 ASSERTL0(false, "fid must be between 0 and 3");
1391 }
1392
1393 nFaceCoeffs = P * (2 * Q - P + 1) / 2;
1394
1395 if (maparray.size() != nFaceCoeffs)
1396 {
1397 maparray = Array<OneD, unsigned int>(nFaceCoeffs);
1398 }
1399
1400 switch (fid)
1401 {
1402 case 0:
1403 idx = 0;
1404 for (i = 0; i < P; ++i)
1405 {
1406 for (j = 0; j < Q - i; ++j)
1407 {
1408 maparray[idx++] = GetMode(i, j, 0);
1409 }
1410 }
1411 break;
1412 case 1:
1413 idx = 0;
1414 for (i = 0; i < P; ++i)
1415 {
1416 for (k = 0; k < Q - i; ++k)
1417 {
1418 maparray[idx++] = GetMode(i, 0, k);
1419 }
1420 }
1421 break;
1422 case 2:
1423 idx = 0;
1424 for (j = 0; j < P - 1; ++j)
1425 {
1426 for (k = 0; k < Q - 1 - j; ++k)
1427 {
1428 maparray[idx++] = GetMode(1, j, k);
1429 // Incorporate modes from zeroth plane where needed.
1430 if (j == 0 && k == 0)
1431 {
1432 maparray[idx++] = GetMode(0, 0, 1);
1433 }
1434 if (j == 0 && k == Q - 2)
1435 {
1436 for (int r = 0; r < Q - 1; ++r)
1437 {
1438 maparray[idx++] = GetMode(0, 1, r);
1439 }
1440 }
1441 }
1442 }
1443 break;
1444 case 3:
1445 idx = 0;
1446 for (j = 0; j < P; ++j)
1447 {
1448 for (k = 0; k < Q - j; ++k)
1449 {
1450 maparray[idx++] = GetMode(0, j, k);
1451 }
1452 }
1453 break;
1454 default:
1455 ASSERTL0(false, "Element map not available.");
1456 }
1457}
1458
1459void StdTetExp::v_GetElmtTraceToTraceMap(const unsigned int fid,
1460 Array<OneD, unsigned int> &maparray,
1461 Array<OneD, int> &signarray,
1462 Orientation faceOrient, int P, int Q)
1463{
1464 int nummodesA = 0, nummodesB = 0, i, j, k, idx;
1465
1467 "Method only implemented for Modified_A BasisType (x "
1468 "direction), Modified_B BasisType (y direction), and "
1469 "Modified_C BasisType(z direction)");
1470
1471 int nFaceCoeffs = 0;
1472
1473 switch (fid)
1474 {
1475 case 0:
1476 nummodesA = m_base[0]->GetNumModes();
1477 nummodesB = m_base[1]->GetNumModes();
1478 break;
1479 case 1:
1480 nummodesA = m_base[0]->GetNumModes();
1481 nummodesB = m_base[2]->GetNumModes();
1482 break;
1483 case 2:
1484 case 3:
1485 nummodesA = m_base[1]->GetNumModes();
1486 nummodesB = m_base[2]->GetNumModes();
1487 break;
1488 default:
1489 ASSERTL0(false, "fid must be between 0 and 3");
1490 }
1491
1492 if (P == -1)
1493 {
1494 P = nummodesA;
1495 Q = nummodesB;
1496 }
1497
1498 nFaceCoeffs = P * (2 * Q - P + 1) / 2;
1499
1500 // Allocate the map array and sign array; set sign array to ones (+)
1501 if (maparray.size() != nFaceCoeffs)
1502 {
1503 maparray = Array<OneD, unsigned int>(nFaceCoeffs, 1);
1504 }
1505
1506 if (signarray.size() != nFaceCoeffs)
1507 {
1508 signarray = Array<OneD, int>(nFaceCoeffs, 1);
1509 }
1510 else
1511 {
1512 fill(signarray.data(), signarray.data() + nFaceCoeffs, 1);
1513 }
1514
1515 // zero signmap and set maparray to zero if elemental
1516 // modes are not as large as face modesl
1517 idx = 0;
1518 int cnt = 0;
1519 int minPA = min(nummodesA, P);
1520 int minQB = min(nummodesB, Q);
1521
1522 for (j = 0; j < minPA; ++j)
1523 {
1524 // set maparray
1525 for (k = 0; k < minQB - j; ++k, ++cnt)
1526 {
1527 maparray[idx++] = cnt;
1528 }
1529
1530 cnt += nummodesB - minQB;
1531
1532 for (k = nummodesB - j; k < Q - j; ++k)
1533 {
1534 signarray[idx] = 0.0;
1535 maparray[idx++] = maparray[0];
1536 }
1537 }
1538
1539 for (j = nummodesA; j < P; ++j)
1540 {
1541 for (k = 0; k < Q - j; ++k)
1542 {
1543 signarray[idx] = 0.0;
1544 maparray[idx++] = maparray[0];
1545 }
1546 }
1547
1548 if (faceOrient == eDir1BwdDir1_Dir2FwdDir2)
1549 {
1550 idx = 0;
1551 for (i = 0; i < P; ++i)
1552 {
1553 for (j = 0; j < Q - i; ++j, idx++)
1554 {
1555 if (i > 1)
1556 {
1557 signarray[idx] = (i % 2 ? -1 : 1);
1558 }
1559 }
1560 }
1561
1562 swap(maparray[0], maparray[Q]);
1563
1564 for (i = 1; i < Q - 1; ++i)
1565 {
1566 swap(maparray[i + 1], maparray[Q + i]);
1567 }
1568 }
1569}
1570
1571/**
1572 * Maps interior modes of an edge to the elemental modes.
1573 */
1575 const int eid, Array<OneD, unsigned int> &maparray,
1576 Array<OneD, int> &signarray, const Orientation edgeOrient)
1577{
1578 int i;
1579 const int P = m_base[0]->GetNumModes();
1580 const int Q = m_base[1]->GetNumModes();
1581 const int R = m_base[2]->GetNumModes();
1582
1583 const int nEdgeIntCoeffs = v_GetEdgeNcoeffs(eid) - 2;
1584
1585 if (maparray.size() != nEdgeIntCoeffs)
1586 {
1587 maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
1588 }
1589 else
1590 {
1591 fill(maparray.data(), maparray.data() + nEdgeIntCoeffs, 0);
1592 }
1593
1594 if (signarray.size() != nEdgeIntCoeffs)
1595 {
1596 signarray = Array<OneD, int>(nEdgeIntCoeffs, 1);
1597 }
1598 else
1599 {
1600 fill(signarray.data(), signarray.data() + nEdgeIntCoeffs, 1);
1601 }
1602
1603 switch (eid)
1604 {
1605 case 0:
1606 for (i = 0; i < P - 2; ++i)
1607 {
1608 maparray[i] = GetMode(i + 2, 0, 0);
1609 }
1610 if (edgeOrient == eBackwards)
1611 {
1612 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1613 {
1614 signarray[i] = -1;
1615 }
1616 }
1617 break;
1618 case 1:
1619 for (i = 0; i < Q - 2; ++i)
1620 {
1621 maparray[i] = GetMode(1, i + 1, 0);
1622 }
1623 if (edgeOrient == eBackwards)
1624 {
1625 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1626 {
1627 signarray[i] = -1;
1628 }
1629 }
1630 break;
1631 case 2:
1632 for (i = 0; i < Q - 2; ++i)
1633 {
1634 maparray[i] = GetMode(0, i + 2, 0);
1635 }
1636 if (edgeOrient == eBackwards)
1637 {
1638 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1639 {
1640 signarray[i] = -1;
1641 }
1642 }
1643 break;
1644 case 3:
1645 for (i = 0; i < R - 2; ++i)
1646 {
1647 maparray[i] = GetMode(0, 0, i + 2);
1648 }
1649 if (edgeOrient == eBackwards)
1650 {
1651 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1652 {
1653 signarray[i] = -1;
1654 }
1655 }
1656 break;
1657 case 4:
1658 for (i = 0; i < R - 2; ++i)
1659 {
1660 maparray[i] = GetMode(1, 0, i + 1);
1661 }
1662 if (edgeOrient == eBackwards)
1663 {
1664 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1665 {
1666 signarray[i] = -1;
1667 }
1668 }
1669 break;
1670 case 5:
1671 for (i = 0; i < R - 2; ++i)
1672 {
1673 maparray[i] = GetMode(0, 1, i + 1);
1674 }
1675 if (edgeOrient == eBackwards)
1676 {
1677 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1678 {
1679 signarray[i] = -1;
1680 }
1681 }
1682 break;
1683 default:
1684 ASSERTL0(false, "Edge not defined.");
1685 break;
1686 }
1687}
1688
1690 const int fid, Array<OneD, unsigned int> &maparray,
1691 Array<OneD, int> &signarray, const Orientation faceOrient)
1692{
1693 int i, j, idx, k;
1694 const int P = m_base[0]->GetNumModes();
1695 const int Q = m_base[1]->GetNumModes();
1696 const int R = m_base[2]->GetNumModes();
1697
1698 const int nFaceIntCoeffs = v_GetTraceIntNcoeffs(fid);
1699
1700 if (maparray.size() != nFaceIntCoeffs)
1701 {
1702 maparray = Array<OneD, unsigned int>(nFaceIntCoeffs);
1703 }
1704
1705 if (signarray.size() != nFaceIntCoeffs)
1706 {
1707 signarray = Array<OneD, int>(nFaceIntCoeffs, 1);
1708 }
1709 else
1710 {
1711 fill(signarray.data(), signarray.data() + nFaceIntCoeffs, 1);
1712 }
1713
1714 switch (fid)
1715 {
1716 case 0:
1717 idx = 0;
1718 for (i = 2; i < P; ++i)
1719 {
1720 for (j = 1; j < Q - i; ++j)
1721 {
1722 if ((int)faceOrient == 7)
1723 {
1724 signarray[idx] = (i % 2 ? -1 : 1);
1725 }
1726 maparray[idx++] = GetMode(i, j, 0);
1727 }
1728 }
1729 break;
1730 case 1:
1731 idx = 0;
1732 for (i = 2; i < P; ++i)
1733 {
1734 for (k = 1; k < R - i; ++k)
1735 {
1736 if ((int)faceOrient == 7)
1737 {
1738 signarray[idx] = (i % 2 ? -1 : 1);
1739 }
1740 maparray[idx++] = GetMode(i, 0, k);
1741 }
1742 }
1743 break;
1744 case 2:
1745 idx = 0;
1746 for (j = 1; j < Q - 1; ++j)
1747 {
1748 for (k = 1; k < R - 1 - j; ++k)
1749 {
1750 if ((int)faceOrient == 7)
1751 {
1752 signarray[idx] = ((j + 1) % 2 ? -1 : 1);
1753 }
1754 maparray[idx++] = GetMode(1, j, k);
1755 }
1756 }
1757 break;
1758 case 3:
1759 idx = 0;
1760 for (j = 2; j < Q; ++j)
1761 {
1762 for (k = 1; k < R - j; ++k)
1763 {
1764 if ((int)faceOrient == 7)
1765 {
1766 signarray[idx] = (j % 2 ? -1 : 1);
1767 }
1768 maparray[idx++] = GetMode(0, j, k);
1769 }
1770 }
1771 break;
1772 default:
1773 ASSERTL0(false, "Face interior map not available.");
1774 break;
1775 }
1776}
1777//---------------------------------------
1778// Wrapper functions
1779//---------------------------------------
1781{
1782
1783 MatrixType mtype = mkey.GetMatrixType();
1784
1785 DNekMatSharedPtr Mat;
1786
1787 switch (mtype)
1788 {
1790 {
1791 int nq0 = m_base[0]->GetNumPoints();
1792 int nq1 = m_base[1]->GetNumPoints();
1793 int nq2 = m_base[2]->GetNumPoints();
1794 int nq;
1795
1796 // take definition from key
1798 {
1799 nq = (int)mkey.GetConstFactor(eFactorConst);
1800 }
1801 else
1802 {
1803 nq = max(nq0, max(nq1, nq2));
1804 }
1805
1806 int neq =
1809 Array<OneD, NekDouble> coll(3);
1811 Array<OneD, NekDouble> tmp(nq0);
1812
1813 Mat =
1814 MemoryManager<DNekMat>::AllocateSharedPtr(neq, nq0 * nq1 * nq2);
1815 int cnt = 0;
1816
1817 for (int i = 0; i < nq; ++i)
1818 {
1819 for (int j = 0; j < nq - i; ++j)
1820 {
1821 for (int k = 0; k < nq - i - j; ++k, ++cnt)
1822 {
1823 coords[cnt] = Array<OneD, NekDouble>(3);
1824 coords[cnt][0] = -1.0 + 2 * k / (NekDouble)(nq - 1);
1825 coords[cnt][1] = -1.0 + 2 * j / (NekDouble)(nq - 1);
1826 coords[cnt][2] = -1.0 + 2 * i / (NekDouble)(nq - 1);
1827 }
1828 }
1829 }
1830
1831 for (int i = 0; i < neq; ++i)
1832 {
1833 LocCoordToLocCollapsed(coords[i], coll);
1834
1835 I[0] = m_base[0]->GetI(coll);
1836 I[1] = m_base[1]->GetI(coll + 1);
1837 I[2] = m_base[2]->GetI(coll + 2);
1838
1839 // interpolate first coordinate direction
1840 NekDouble fac;
1841 for (int k = 0; k < nq2; ++k)
1842 {
1843 for (int j = 0; j < nq1; ++j)
1844 {
1845
1846 fac = (I[1]->GetPtr())[j] * (I[2]->GetPtr())[k];
1847 Vmath::Smul(nq0, fac, I[0]->GetPtr(), 1, tmp, 1);
1848
1849 Vmath::Vcopy(nq0, &tmp[0], 1,
1850 Mat->GetRawPtr() + k * nq0 * nq1 * neq +
1851 j * nq0 * neq + i,
1852 neq);
1853 }
1854 }
1855 }
1856 }
1857 break;
1858 default:
1859 {
1861 }
1862 break;
1863 }
1864
1865 return Mat;
1866}
1867
1869{
1870 return v_GenMatrix(mkey);
1871}
1872
1873//---------------------------------------
1874// Private helper functions
1875//---------------------------------------
1876
1877/**
1878 * @brief Compute the mode number in the expansion for a particular
1879 * tensorial combination.
1880 *
1881 * Modes are numbered with the r index travelling fastest, followed by
1882 * q and then p, and each q-r plane is of size
1883 * (Q+1)*(Q+2)/2+max(0,R-Q-p)*Q. For example, when P=2, Q=3 and R=4 (nm0=3, nm1
1884 * = 4, nm2 = 5) the indexing inside each q-r plane (with r increasing upwards
1885 * and q to the right) is:
1886 *
1887 * 4
1888 * 3 8 17
1889 * 2 7 11 16 20 25
1890 * 1 6 10 13 15 19 22 24 27
1891 * 0 5 9 12 14 18 21 23 26
1892 *
1893 * Geometrically they can be interpreted as
1894 * p = 0: p = 2: p = 1:
1895 * ----------------------------------
1896 * 1
1897 * 4 8 17
1898 * 3 11 7 25 16 20
1899 * 2 10 13 6 24 27 15 19 22
1900 * 0 9 12 5 23 26 14 18 21
1901 *
1902 * so we have the following breakdown
1903 *
1904 * Vertices V[0,1,2,3] = [0, 14, 5, 1]
1905 * Edges E[0,1,2,3,4,5,6] =[[23],[18, 21],
1906 * [9, 12], [2,3,4], [15, 16, 17], [6,7,8]]
1907 * Faces F[0.1,2,3] = [[26], [24,25], [19, 22, 20], [10, 13, 11]
1908 * Interior [27]
1909 * Note that in this element, we must have that \f$ P \leq Q \leq
1910 * R\f$.
1911 */
1912int StdTetExp::GetMode(const int I, const int J, const int K)
1913{
1914 const int Q = m_base[1]->GetNumModes();
1915 const int R = m_base[2]->GetNumModes();
1916
1917 int i, j, q_hat, k_hat;
1918 int cnt = 0;
1919
1920 // Traverse to q-r plane number I
1921 for (i = 0; i < I; ++i)
1922 {
1923 // Size of triangle part
1924 q_hat = Q - i;
1925 // Size of rectangle part
1926 k_hat = R - Q;
1927 cnt += q_hat * (q_hat + 1) / 2 + k_hat * (Q - i);
1928 }
1929
1930 // Traverse to q column J
1931 q_hat = R - I;
1932 for (j = 0; j < J; ++j)
1933 {
1934 cnt += q_hat;
1935 q_hat--;
1936 }
1937
1938 // Traverse up stacks to K
1939 cnt += K;
1940
1941 return cnt;
1942}
1943
1945 const Array<OneD, const NekDouble> &inarray,
1946 Array<OneD, NekDouble> &outarray)
1947{
1948 int i, j;
1949
1950 int nquad0 = m_base[0]->GetNumPoints();
1951 int nquad1 = m_base[1]->GetNumPoints();
1952 int nquad2 = m_base[2]->GetNumPoints();
1953
1954 const Array<OneD, const NekDouble> &w0 = m_base[0]->GetW();
1955 const Array<OneD, const NekDouble> &w1 = m_base[1]->GetW();
1956 const Array<OneD, const NekDouble> &w2 = m_base[2]->GetW();
1957
1958 const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
1959 const Array<OneD, const NekDouble> &z2 = m_base[2]->GetZ();
1960
1961 // multiply by integration constants
1962 for (i = 0; i < nquad1 * nquad2; ++i)
1963 {
1964 Vmath::Vmul(nquad0, (NekDouble *)&inarray[0] + i * nquad0, 1, w0.data(),
1965 1, &outarray[0] + i * nquad0, 1);
1966 }
1967
1968 switch (m_base[1]->GetPointsType())
1969 {
1970 // (1,0) Jacobi Inner product.
1971 case LibUtilities::eGaussRadauMAlpha1Beta0:
1972 for (j = 0; j < nquad2; ++j)
1973 {
1974 for (i = 0; i < nquad1; ++i)
1975 {
1976 Blas::Dscal(nquad0, 0.5 * w1[i],
1977 &outarray[0] + i * nquad0 + j * nquad0 * nquad1,
1978 1);
1979 }
1980 }
1981 break;
1982
1983 default:
1984 for (j = 0; j < nquad2; ++j)
1985 {
1986 for (i = 0; i < nquad1; ++i)
1987 {
1988 Blas::Dscal(nquad0, 0.5 * (1 - z1[i]) * w1[i],
1989 &outarray[0] + i * nquad0 + j * nquad0 * nquad1,
1990 1);
1991 }
1992 }
1993 break;
1994 }
1995
1996 switch (m_base[2]->GetPointsType())
1997 {
1998 // (2,0) Jacobi inner product.
1999 case LibUtilities::eGaussRadauMAlpha2Beta0:
2000 for (i = 0; i < nquad2; ++i)
2001 {
2002 Blas::Dscal(nquad0 * nquad1, 0.25 * w2[i],
2003 &outarray[0] + i * nquad0 * nquad1, 1);
2004 }
2005 break;
2006 // (1,0) Jacobi inner product.
2007 case LibUtilities::eGaussRadauMAlpha1Beta0:
2008 for (i = 0; i < nquad2; ++i)
2009 {
2010 Blas::Dscal(nquad0 * nquad1, 0.25 * (1 - z2[i]) * w2[i],
2011 &outarray[0] + i * nquad0 * nquad1, 1);
2012 }
2013 break;
2014 default:
2015 for (i = 0; i < nquad2; ++i)
2016 {
2017 Blas::Dscal(nquad0 * nquad1,
2018 0.25 * (1 - z2[i]) * (1 - z2[i]) * w2[i],
2019 &outarray[0] + i * nquad0 * nquad1, 1);
2020 }
2021 break;
2022 }
2023}
2024
2026 const StdMatrixKey &mkey)
2027{
2028 // To do : 1) add a test to ensure 0 \leq SvvCutoff \leq 1.
2029 // 2) check if the transfer function needs an analytical
2030 // Fourier transform.
2031 // 3) if it doesn't : find a transfer function that renders
2032 // the if( cutoff_a ...) useless to reduce computational
2033 // cost.
2034 // 4) add SVVDiffCoef to both models!!
2035
2036 int qa = m_base[0]->GetNumPoints();
2037 int qb = m_base[1]->GetNumPoints();
2038 int qc = m_base[2]->GetNumPoints();
2039 int nmodes_a = m_base[0]->GetNumModes();
2040 int nmodes_b = m_base[1]->GetNumModes();
2041 int nmodes_c = m_base[2]->GetNumModes();
2042
2043 // Declare orthogonal basis.
2047
2051
2052 StdTetExp OrthoExp(Ba, Bb, Bc);
2053
2054 Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
2055 int i, j, k, cnt = 0;
2056
2057 // project onto physical space.
2058 OrthoExp.FwdTrans(array, orthocoeffs);
2059
2061 {
2062 // Rodrigo's power kernel
2064 NekDouble SvvDiffCoeff =
2067
2068 for (i = 0; i < nmodes_a; ++i)
2069 {
2070 for (j = 0; j < nmodes_b - j; ++j)
2071 {
2072 NekDouble fac1 = std::max(
2073 pow((1.0 * i) / (nmodes_a - 1), cutoff * nmodes_a),
2074 pow((1.0 * j) / (nmodes_b - 1), cutoff * nmodes_b));
2075
2076 for (k = 0; k < nmodes_c - i - j; ++k)
2077 {
2078 NekDouble fac =
2079 std::max(fac1, pow((1.0 * k) / (nmodes_c - 1),
2080 cutoff * nmodes_c));
2081
2082 orthocoeffs[cnt] *= SvvDiffCoeff * fac;
2083 cnt++;
2084 }
2085 }
2086 }
2087 }
2088 else if (mkey.ConstFactorExists(
2089 eFactorSVVDGKerDiffCoeff)) // Rodrigo/Mansoor's DG Kernel
2090 {
2093
2094 int max_abc = max(nmodes_a - kSVVDGFiltermodesmin,
2095 nmodes_b - kSVVDGFiltermodesmin);
2096 max_abc = max(max_abc, nmodes_c - kSVVDGFiltermodesmin);
2097 // clamp max_abc
2098 max_abc = max(max_abc, 0);
2099 max_abc = min(max_abc, kSVVDGFiltermodesmax - kSVVDGFiltermodesmin);
2100
2101 for (i = 0; i < nmodes_a; ++i)
2102 {
2103 for (j = 0; j < nmodes_b - j; ++j)
2104 {
2105 int maxij = max(i, j);
2106
2107 for (k = 0; k < nmodes_c - i - j; ++k)
2108 {
2109 int maxijk = max(maxij, k);
2110 maxijk = min(maxijk, kSVVDGFiltermodesmax - 1);
2111
2112 orthocoeffs[cnt] *=
2113 SvvDiffCoeff * kSVVDGFilter[max_abc][maxijk];
2114 cnt++;
2115 }
2116 }
2117 }
2118 }
2119 else
2120 {
2121
2122 // SVV filter paramaters (how much added diffusion
2123 // relative to physical one and fraction of modes from
2124 // which you start applying this added diffusion)
2125
2126 NekDouble SvvDiffCoeff =
2128 NekDouble SVVCutOff =
2130
2131 // Defining the cut of mode
2132 int cutoff_a = (int)(SVVCutOff * nmodes_a);
2133 int cutoff_b = (int)(SVVCutOff * nmodes_b);
2134 int cutoff_c = (int)(SVVCutOff * nmodes_c);
2135 int nmodes = min(min(nmodes_a, nmodes_b), nmodes_c);
2136 NekDouble cutoff = min(min(cutoff_a, cutoff_b), cutoff_c);
2137 NekDouble epsilon = 1;
2138
2139 //------"New" Version August 22nd '13--------------------
2140 for (i = 0; i < nmodes_a; ++i)
2141 {
2142 for (j = 0; j < nmodes_b - i; ++j)
2143 {
2144 for (k = 0; k < nmodes_c - i - j; ++k)
2145 {
2146 if (i + j + k >= cutoff)
2147 {
2148 orthocoeffs[cnt] *= ((SvvDiffCoeff)*exp(
2149 -(i + j + k - nmodes) * (i + j + k - nmodes) /
2150 ((NekDouble)((i + j + k - cutoff + epsilon) *
2151 (i + j + k - cutoff + epsilon)))));
2152 }
2153 else
2154 {
2155 orthocoeffs[cnt] *= 0.0;
2156 }
2157 cnt++;
2158 }
2159 }
2160 }
2161 }
2162
2163 // backward transform to physical space
2164 OrthoExp.BwdTrans(orthocoeffs, array);
2165}
2166
2168 const Array<OneD, const NekDouble> &inarray,
2169 Array<OneD, NekDouble> &outarray)
2170{
2171 int nquad0 = m_base[0]->GetNumPoints();
2172 int nquad1 = m_base[1]->GetNumPoints();
2173 int nquad2 = m_base[2]->GetNumPoints();
2174 int nqtot = nquad0 * nquad1 * nquad2;
2175 int nmodes0 = m_base[0]->GetNumModes();
2176 int nmodes1 = m_base[1]->GetNumModes();
2177 int nmodes2 = m_base[2]->GetNumModes();
2178 int numMax = nmodes0;
2179
2181 Array<OneD, NekDouble> coeff_tmp1(m_ncoeffs, 0.0);
2182 Array<OneD, NekDouble> coeff_tmp2(m_ncoeffs, 0.0);
2183 Array<OneD, NekDouble> phys_tmp(nqtot, 0.0);
2184 Array<OneD, NekDouble> tmp, tmp2, tmp3, tmp4;
2185
2186 Vmath::Vcopy(m_ncoeffs, inarray, 1, coeff_tmp2, 1);
2187
2188 const LibUtilities::PointsKey Pkey0 = m_base[0]->GetPointsKey();
2189 const LibUtilities::PointsKey Pkey1 = m_base[1]->GetPointsKey();
2190 const LibUtilities::PointsKey Pkey2 = m_base[2]->GetPointsKey();
2191
2192 LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A, nmodes0, Pkey0);
2193 LibUtilities::BasisKey bortho1(LibUtilities::eOrtho_B, nmodes1, Pkey1);
2194 LibUtilities::BasisKey bortho2(LibUtilities::eOrtho_C, nmodes2, Pkey2);
2195
2196 Vmath::Zero(m_ncoeffs, coeff_tmp2, 1);
2197
2200 bortho0, bortho1, bortho2);
2201
2202 BwdTrans(inarray, phys_tmp);
2203 OrthoTetExp->FwdTrans(phys_tmp, coeff);
2204
2205 Vmath::Zero(m_ncoeffs, outarray, 1);
2206
2207 // filtering
2208 int cnt = 0;
2209 for (int u = 0; u < numMin; ++u)
2210 {
2211 for (int i = 0; i < numMin - u; ++i)
2212 {
2213 Vmath::Vcopy(numMin - u - i, tmp = coeff + cnt, 1,
2214 tmp2 = coeff_tmp1 + cnt, 1);
2215 cnt += numMax - u - i;
2216 }
2217 for (int i = numMin; i < numMax - u; ++i)
2218 {
2219 cnt += numMax - u - i;
2220 }
2221 }
2222
2223 OrthoTetExp->BwdTrans(coeff_tmp1, phys_tmp);
2224 FwdTrans(phys_tmp, outarray);
2225}
2226
2228 Array<OneD, int> &conn, [[maybe_unused]] bool standard)
2229{
2230 int np0 = m_base[0]->GetNumPoints();
2231 int np1 = m_base[1]->GetNumPoints();
2232 int np2 = m_base[2]->GetNumPoints();
2233 int np = max(np0, max(np1, np2));
2234
2235 conn = Array<OneD, int>(4 * (np - 1) * (np - 1) * (np - 1));
2236
2237 int row = 0;
2238 int rowp1 = 0;
2239 int plane = 0;
2240 int row1 = 0;
2241 int row1p1 = 0;
2242 int planep1 = 0;
2243 int cnt = 0;
2244 for (int i = 0; i < np - 1; ++i)
2245 {
2246 planep1 += (np - i) * (np - i + 1) / 2;
2247 row = 0; // current plane row offset
2248 rowp1 = 0; // current plane row plus one offset
2249 row1 = 0; // next plane row offset
2250 row1p1 = 0; // nex plane row plus one offset
2251 for (int j = 0; j < np - i - 1; ++j)
2252 {
2253 rowp1 += np - i - j;
2254 row1p1 += np - i - j - 1;
2255 for (int k = 0; k < np - i - j - 2; ++k)
2256 {
2257 conn[cnt++] = plane + row + k + 1;
2258 conn[cnt++] = plane + row + k;
2259 conn[cnt++] = plane + rowp1 + k;
2260 conn[cnt++] = planep1 + row1 + k;
2261
2262 conn[cnt++] = plane + row + k + 1;
2263 conn[cnt++] = plane + rowp1 + k + 1;
2264 conn[cnt++] = planep1 + row1 + k + 1;
2265 conn[cnt++] = planep1 + row1 + k;
2266
2267 conn[cnt++] = plane + rowp1 + k + 1;
2268 conn[cnt++] = plane + row + k + 1;
2269 conn[cnt++] = plane + rowp1 + k;
2270 conn[cnt++] = planep1 + row1 + k;
2271
2272 conn[cnt++] = planep1 + row1 + k;
2273 conn[cnt++] = planep1 + row1p1 + k;
2274 conn[cnt++] = plane + rowp1 + k;
2275 conn[cnt++] = plane + rowp1 + k + 1;
2276
2277 conn[cnt++] = planep1 + row1 + k;
2278 conn[cnt++] = planep1 + row1p1 + k;
2279 conn[cnt++] = planep1 + row1 + k + 1;
2280 conn[cnt++] = plane + rowp1 + k + 1;
2281
2282 if (k < np - i - j - 3)
2283 {
2284 conn[cnt++] = plane + rowp1 + k + 1;
2285 conn[cnt++] = planep1 + row1p1 + k + 1;
2286 conn[cnt++] = planep1 + row1 + k + 1;
2287 conn[cnt++] = planep1 + row1p1 + k;
2288 }
2289 }
2290
2291 conn[cnt++] = plane + row + np - i - j - 1;
2292 conn[cnt++] = plane + row + np - i - j - 2;
2293 conn[cnt++] = plane + rowp1 + np - i - j - 2;
2294 conn[cnt++] = planep1 + row1 + np - i - j - 2;
2295
2296 row += np - i - j;
2297 row1 += np - i - j - 1;
2298 }
2299 plane += (np - i) * (np - i + 1) / 2;
2300 }
2301}
2302
2303} // 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
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 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.
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.
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Backward transformation from coefficient space to physical space.
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.
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
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
MatrixType GetMatrixType() const
NekDouble GetConstFactor(const ConstFactorType &factor) const
bool ConstFactorExists(const ConstFactorType &factor) const
int v_GetNtraces() const override
int v_GetTraceNcoeffs(const int i) const override
NekDouble v_PhysEvalFirstDeriv(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
void v_GetTraceNumModes(const int fid, int &numModes0, int &numModes1, Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2) override
int v_GetTraceIntNcoeffs(const int i) const override
void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray) override
void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_GetInteriorMap(Array< OneD, unsigned int > &outarray) override
int v_GetNedges() const override
int GetMode(const int i, const int j, const int k)
Compute the mode number in the expansion for a particular tensorial combination.
int v_GetEdgeNcoeffs(const int i) const override
bool v_IsBoundaryInteriorExpansion() const override
LibUtilities::PointsKey v_GetTracePointsKey(const int i, const int j) const override
NekDouble v_PhysEvaluateBasis(const Array< OneD, const NekDouble > &coords, int mode) final
void v_GetEdgeInteriorToElementMap(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
StdTetExp(const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
Definition StdTetExp.cpp:42
int v_GetTraceNumPoints(const int i) const override
void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) 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
void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey) override
DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey) override
void v_LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta) override
void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray) override
void v_GetSimplexEquiSpacedConnectivity(Array< OneD, int > &conn, bool standard=true) override
void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
int v_NumBndryCoeffs() const override
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) 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
void v_GetCoords(Array< OneD, NekDouble > &coords_x, Array< OneD, NekDouble > &coords_y, Array< OneD, NekDouble > &coords_z) override
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_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_GetElmtTraceToTraceMap(const unsigned int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1) override
int v_GetVertexMap(int localVertexId, bool useCoeffPacking=false) override
int v_GetNverts() const override
void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dx, Array< OneD, NekDouble > &out_dy, Array< OneD, NekDouble > &out_dz) override
Calculate the derivative of the physical points.
Definition StdTetExp.cpp:83
const LibUtilities::BasisKey v_GetTraceBasisKey(const int i, const int k, bool UseGLL=false) const override
LibUtilities::ShapeType v_DetShapeType() const override
void v_LocCollapsedToLocCoord(const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi) override
LibUtilities::ShapeType DetShapeType() const
Definition StdTetExp.h:56
int v_NumDGBndryCoeffs() const 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
int v_CalcNumberOfCoefficients(const std::vector< unsigned int > &nummodes, int &modes_offset) override
DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) 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 getNumberOfBndCoefficients(int Na, int Nb, int Nc)
int getNumberOfCoefficients(int Na, int Nb, int Nc)
@ eModified_B
Principle Modified Functions .
Definition BasisType.h:49
@ P
Monomial polynomials .
Definition BasisType.h:62
@ eOrtho_A
Principle Orthogonal Functions .
Definition BasisType.h:42
@ eModified_C
Principle Modified Functions .
Definition BasisType.h:50
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition BasisType.h:56
@ eOrtho_C
Principle Orthogonal Functions .
Definition BasisType.h:46
@ eOrtho_B
Principle Orthogonal Functions .
Definition BasisType.h:44
@ eModified_A
Principle Modified Functions .
Definition BasisType.h:48
static const NekDouble kNekZeroTol
std::shared_ptr< StdTetExp > StdTetExpSharedPtr
Definition StdTetExp.h:224
LibUtilities::BasisKey EvaluateTriFaceBasisKey(const int facedir, const LibUtilities::BasisSharedPtr &faceDirBasis, bool UseGLL)
const int kSVVDGFiltermodesmin
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 Zero(int n, T *x, const int incx)
Zero vector.
Definition Vmath.hpp:273
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition Vmath.hpp:194
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
scalarT< T > sqrt(scalarT< T > in)
Definition scalar.hpp:290