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