47 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
50 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
55 "order in 'a' direction is higher "
56 "than order in 'c' direction");
58 "order in 'b' direction is higher "
59 "than order in 'c' direction");
62 "Expected basis type in 'c' direction to be ModifiedPyr_C or "
89 int Qx =
m_base[0]->GetNumPoints();
90 int Qy =
m_base[1]->GetNumPoints();
91 int Qz =
m_base[2]->GetNumPoints();
100 eta_y =
m_base[1]->GetZ();
101 eta_z =
m_base[2]->GetZ();
105 if (out_dxi1.size() > 0)
107 for (k = 0, n = 0; k < Qz; ++k)
110 for (j = 0; j < Qy; ++j)
112 for (i = 0; i < Qx; ++i, ++n)
114 out_dxi1[n] = fac * dEta_bar1[n];
120 if (out_dxi2.size() > 0)
122 for (k = 0, n = 0; k < Qz; ++k)
125 for (j = 0; j < Qy; ++j)
127 for (i = 0; i < Qx; ++i, ++n)
129 out_dxi2[n] = fac * dXi2[n];
135 if (out_dxi3.size() > 0)
137 for (k = 0, n = 0; k < Qz; ++k)
140 for (j = 0; j < Qy; ++j)
143 for (i = 0; i < Qx; ++i, ++n)
145 out_dxi3[n] = (1.0 + eta_x[i]) * fac * dEta_bar1[n] +
146 fac1 * fac * dXi2[n] + dEta3[n];
182 ASSERTL1(
false,
"input dir is out of range");
231 if (
m_base[0]->Collocation() &&
m_base[1]->Collocation() &&
236 inarray, 1, outarray, 1);
250 int nquad0 =
m_base[0]->GetNumPoints();
251 int nquad1 =
m_base[1]->GetNumPoints();
252 int nquad2 =
m_base[2]->GetNumPoints();
253 int order0 =
m_base[0]->GetNumModes();
254 int order1 =
m_base[1]->GetNumModes();
257 nquad2 * nquad1 * nquad0);
260 m_base[2]->GetBdata(), inarray, outarray, wsp,
true,
270 [[maybe_unused]]
bool doCheckCollDir0,
271 [[maybe_unused]]
bool doCheckCollDir1,
272 [[maybe_unused]]
bool doCheckCollDir2)
274 int nquad0 =
m_base[0]->GetNumPoints();
275 int nquad1 =
m_base[1]->GetNumPoints();
276 int nquad2 =
m_base[2]->GetNumPoints();
278 int order0 =
m_base[0]->GetNumModes();
279 int order1 =
m_base[1]->GetNumModes();
280 int order2 =
m_base[2]->GetNumModes();
285 int i, j, mode, mode1, cnt;
288 mode = mode1 = cnt = 0;
289 for (i = 0; i < order0; ++i)
291 for (j = 0; j < order1; ++j, ++cnt)
293 int ijmax = max(i, j);
295 base2.get() + mode * nquad2, nquad2,
296 inarray.get() + mode1, 1, 0.0, tmp.get() + cnt * nquad2,
298 mode += order2 - ijmax;
299 mode1 += order2 - ijmax;
302 for (j = order1; j < order2; ++j)
304 int ijmax = max(i, j);
305 mode += order2 - ijmax;
317 Blas::Daxpy(nquad2, inarray[1], base2.get() + nquad2, 1,
318 &tmp[0] + nquad2, 1);
321 Blas::Daxpy(nquad2, inarray[1], base2.get() + nquad2, 1,
322 &tmp[0] + order1 * nquad2, 1);
325 Blas::Daxpy(nquad2, inarray[1], base2.get() + nquad2, 1,
326 &tmp[0] + order1 * nquad2 + nquad2, 1);
331 for (i = 0; i < order0; ++i)
333 Blas::Dgemm(
'N',
'T', nquad1, nquad2, order1, 1.0, base1.get(), nquad1,
334 tmp.get() + mode * nquad2, nquad2, 0.0,
335 tmp1.get() + i * nquad1 * nquad2, nquad1);
340 Blas::Dgemm(
'N',
'T', nquad0, nquad1 * nquad2, order0, 1.0, base0.get(),
341 nquad0, tmp1.get(), nquad1 * nquad2, 0.0, outarray.get(),
371 out = (*imatsys) * in;
396 if (
m_base[0]->Collocation() &&
m_base[1]->Collocation() &&
412 int nquad1 =
m_base[1]->GetNumPoints();
413 int nquad2 =
m_base[2]->GetNumPoints();
414 int order0 =
m_base[0]->GetNumModes();
415 int order1 =
m_base[1]->GetNumModes();
419 if (multiplybyweights)
427 tmp, outarray, wsp,
true,
true,
true);
433 inarray, outarray, wsp,
true,
true,
true);
443 [[maybe_unused]]
bool doCheckCollDir0,
444 [[maybe_unused]]
bool doCheckCollDir1,
445 [[maybe_unused]]
bool doCheckCollDir2)
447 int nquad0 =
m_base[0]->GetNumPoints();
448 int nquad1 =
m_base[1]->GetNumPoints();
449 int nquad2 =
m_base[2]->GetNumPoints();
451 int order0 =
m_base[0]->GetNumModes();
452 int order1 =
m_base[1]->GetNumModes();
453 int order2 =
m_base[2]->GetNumModes();
455 ASSERTL1(wsp.size() >= nquad1 * nquad2 * order0 + nquad2 * order0 * order1,
456 "Insufficient workspace size");
461 int i, j, mode, mode1, cnt;
464 Blas::Dgemm(
'T',
'N', nquad1 * nquad2, order0, nquad0, 1.0, inarray.get(),
465 nquad0, base0.get(), nquad0, 0.0, tmp1.get(), nquad1 * nquad2);
468 for (mode = i = 0; i < order0; ++i)
471 tmp1.get() + i * nquad1 * nquad2, nquad1, base1.get(),
472 nquad1, 0.0, tmp2.get() + mode * nquad2, nquad2);
477 mode = mode1 = cnt = 0;
478 for (i = 0; i < order0; ++i)
480 for (j = 0; j < order1; ++j, ++cnt)
482 int ijmax = max(i, j);
485 base2.get() + mode * nquad2, nquad2,
486 tmp2.get() + cnt * nquad2, 1, 0.0,
487 outarray.get() + mode1, 1);
488 mode += order2 - ijmax;
489 mode1 += order2 - ijmax;
493 for (j = order1; j < order2; ++j)
495 int ijmax = max(i, j);
496 mode += order2 - ijmax;
506 Blas::Ddot(nquad2, base2.get() + nquad2, 1, &tmp2[nquad2], 1);
509 outarray[1] +=
Blas::Ddot(nquad2, base2.get() + nquad2, 1,
510 &tmp2[nquad2 * order1], 1);
513 outarray[1] +=
Blas::Ddot(nquad2, base2.get() + nquad2, 1,
514 &tmp2[nquad2 * order1 + nquad2], 1);
536 int nquad0 =
m_base[0]->GetNumPoints();
537 int nquad1 =
m_base[1]->GetNumPoints();
538 int nquad2 =
m_base[2]->GetNumPoints();
539 int nqtot = nquad0 * nquad1 * nquad2;
546 int order0 =
m_base[0]->GetNumModes();
547 int order1 =
m_base[1]->GetNumModes();
550 nquad2 * order0 * order1);
557 for (i = 0; i < nquad0; ++i)
559 gfac0[i] = 0.5 * (1 + z0[i]);
563 for (i = 0; i < nquad1; ++i)
565 gfac1[i] = 0.5 * (1 + z1[i]);
569 for (i = 0; i < nquad2; ++i)
571 gfac2[i] = 2.0 / (1 - z2[i]);
575 const int nq01 = nquad0 * nquad1;
576 for (i = 0; i < nquad2; ++i)
578 Vmath::Smul(nq01, gfac2[i], &inarray[0] + i * nq01, 1,
579 &tmp0[0] + i * nq01, 1);
590 m_base[2]->GetBdata(), tmp0, outarray, wsp,
false,
true,
true);
597 m_base[2]->GetBdata(), tmp0, outarray, wsp,
true,
false,
true);
606 for (i = 0; i < nquad1 * nquad2; ++i)
608 Vmath::Vmul(nquad0, tmp0.get() + i * nquad0, 1, gfac0.get(), 1,
609 tmp0.get() + i * nquad0, 1);
613 m_base[2]->GetBdata(), tmp0, tmp3, wsp,
false,
true,
true);
616 for (i = 0; i < nquad2; ++i)
618 Vmath::Smul(nq01, gfac2[i], &inarray[0] + i * nq01, 1,
619 &tmp0[0] + i * nq01, 1);
621 for (i = 0; i < nquad1 * nquad2; ++i)
623 Vmath::Smul(nquad0, gfac1[i % nquad1], &tmp0[0] + i * nquad0, 1,
624 &tmp0[0] + i * nquad0, 1);
630 m_base[2]->GetBdata(), tmp0, tmp4, wsp,
true,
false,
true);
635 m_base[2]->GetDbdata(), tmp0, outarray, wsp,
true,
true,
false);
645 ASSERTL1(
false,
"input dir is out of range");
671 eta[1] = 2.0 * (1.0 + xi[1]) / d2 - 1.0;
672 eta[0] = 2.0 * (1.0 + xi[0]) / d2 - 1.0;
678 xi[0] = (1.0 + eta[0]) * (1.0 - eta[2]) * 0.5 - 1.0;
679 xi[1] = (1.0 + eta[1]) * (1.0 - eta[2]) * 0.5 - 1.0;
695 for (
int k = 0; k < Qz; ++k)
697 for (
int j = 0; j < Qy; ++j)
699 for (
int i = 0; i < Qx; ++i)
701 int s = i + Qx * (j + Qy * k);
704 xi_y[s] = (1.0 + eta_y[j]) * (1.0 - eta_z[k]) / 2.0 - 1.0;
705 xi_x[s] = (1.0 + etaBar_x[i]) * (1.0 - eta_z[k]) / 2.0 - 1.0;
713 std::array<NekDouble, 3> &firstOrderDerivs)
720 if ((1 - coll[2]) < 1e-5)
724 EphysDeriv2(totPoints);
725 PhysDeriv(inarray, EphysDeriv0, EphysDeriv1, EphysDeriv2);
728 I[0] =
GetBase()[0]->GetI(coll);
729 I[1] =
GetBase()[1]->GetI(coll + 1);
730 I[2] =
GetBase()[2]->GetI(coll + 2);
738 std::array<NekDouble, 3> interDeriv;
743 firstOrderDerivs[0] = fac * interDeriv[0];
744 firstOrderDerivs[1] = fac * interDeriv[1];
745 firstOrderDerivs[2] = ((1.0 + coll[0]) / (1.0 - coll[2])) * interDeriv[0] +
746 ((1.0 + coll[1]) / (1.0 - coll[2])) * interDeriv[1] +
762 int nummodes[3] = {
m_base[0]->GetNumModes(),
m_base[1]->GetNumModes(),
763 m_base[2]->GetNumModes()};
769 numModes0 = nummodes[0];
770 numModes1 = nummodes[1];
776 numModes0 = nummodes[0];
777 numModes1 = nummodes[2];
783 numModes0 = nummodes[1];
784 numModes1 = nummodes[2];
791 std::swap(numModes0, numModes1);
801 const int nm0 =
m_base[0]->GetNumModes();
802 const int nm1 =
m_base[1]->GetNumModes();
803 const int nm2 =
m_base[2]->GetNumModes();
805 int mode0 = 0, mode1 = 0, mode2 = 0, cnt = 0;
808 for (mode0 = 0; mode0 < nm0; ++mode0)
810 for (mode1 = 0; mode1 < nm1; ++mode1)
812 int maxpq = max(mode0, mode1);
813 for (mode2 = 0; mode2 < nm2 - maxpq; ++mode2, ++cnt)
833 for (
int j = nm1; j < nm2; ++j)
835 int ijmax = max(mode0, j);
836 mode2 += nm2 - ijmax;
842 return StdExpansion::BaryEvaluateBasis<2>(coll[2], 1);
846 return StdExpansion::BaryEvaluateBasis<0>(coll[0], mode0) *
847 StdExpansion::BaryEvaluateBasis<1>(coll[1], mode1) *
848 StdExpansion::BaryEvaluateBasis<2>(coll[2], mode2);
880 "BasisType is not a boundary interior form");
883 "BasisType is not a boundary interior form");
886 "BasisType is not a boundary interior form");
888 int P =
m_base[0]->GetNumModes();
889 int Q =
m_base[1]->GetNumModes();
890 int R =
m_base[2]->GetNumModes();
899 "BasisType is not a boundary interior form");
902 "BasisType is not a boundary interior form");
905 "BasisType is not a boundary interior form");
907 int P =
m_base[0]->GetNumModes() - 1;
908 int Q =
m_base[1]->GetNumModes() - 1;
909 int R =
m_base[2]->GetNumModes() - 1;
911 return (
P + 1) * (Q + 1)
912 + 2 * (R + 1) +
P * (1 + 2 * R -
P)
913 + 2 * (R + 1) + Q * (1 + 2 * R - Q);
918 ASSERTL2(i >= 0 && i <= 4,
"face id is out of range");
924 else if (i == 1 || i == 3)
927 return Q + 1 + (
P * (1 + 2 * Q -
P)) / 2;
932 return Q + 1 + (
P * (1 + 2 * Q -
P)) / 2;
938 ASSERTL2(i >= 0 && i <= 4,
"face id is out of range");
940 int P =
m_base[0]->GetNumModes() - 1;
941 int Q =
m_base[1]->GetNumModes() - 1;
942 int R =
m_base[2]->GetNumModes() - 1;
946 return (
P - 1) * (Q - 1);
948 else if (i == 1 || i == 3)
950 return (
P - 1) * (2 * (R - 1) - (
P - 1) - 1) / 2;
954 return (Q - 1) * (2 * (R - 1) - (Q - 1) - 1) / 2;
960 ASSERTL2(i >= 0 && i <= 4,
"face id is out of range");
964 return m_base[0]->GetNumPoints() *
m_base[1]->GetNumPoints();
966 else if (i == 1 || i == 3)
968 return m_base[0]->GetNumPoints() *
m_base[2]->GetNumPoints();
972 return m_base[1]->GetNumPoints() *
m_base[2]->GetNumPoints();
978 ASSERTL2(i >= 0 && i <= 7,
"edge id is out of range");
980 if (i == 0 || i == 2)
984 else if (i == 1 || i == 3)
997 ASSERTL2(i >= 0 && i <= 4,
"face id is out of range");
998 ASSERTL2(k >= 0 && k <= 1,
"basis key id is out of range");
1006 m_base[k]->GetNumModes());
1013 m_base[2 * k]->GetNumModes());
1020 m_base[k + 1]->GetNumModes());
1029 const std::vector<unsigned int> &nummodes,
int &modes_offset)
1032 nummodes[modes_offset], nummodes[modes_offset + 1],
1033 nummodes[modes_offset + 2]);
1044 "Mapping not defined for this type of basis");
1048 if (useCoeffPacking ==
true)
1068 ASSERTL0(
false,
"local vertex id must be between 0 and 4");
1091 ASSERTL0(
false,
"local vertex id must be between 0 and 4");
1102 "BasisType is not a boundary interior form");
1105 "BasisType is not a boundary interior form");
1108 "BasisType is not a boundary interior form");
1110 int P =
m_base[0]->GetNumModes() - 1,
p;
1111 int Q =
m_base[1]->GetNumModes() - 1,
q;
1112 int R =
m_base[2]->GetNumModes() - 1, r;
1116 if (outarray.size() != nIntCoeffs)
1124 for (
p = 2;
p <=
P; ++
p)
1126 for (
q = 2;
q <= Q; ++
q)
1128 int maxpq = max(
p,
q);
1129 for (r = 1; r <= R - maxpq; ++r)
1141 "BasisType is not a boundary interior form");
1144 "BasisType is not a boundary interior form");
1147 "BasisType is not a boundary interior form");
1149 int P =
m_base[0]->GetNumModes() - 1,
p;
1150 int Q =
m_base[1]->GetNumModes() - 1,
q;
1151 int R =
m_base[2]->GetNumModes() - 1, r;
1156 if (maparray.size() != nBnd)
1162 for (
p = 0;
p <=
P; ++
p)
1167 for (
q = 0;
q <= Q; ++
q)
1169 int maxpq = max(
p,
q);
1170 for (r = 0; r <= R - maxpq; ++r)
1180 for (
q = 0;
q <= Q; ++
q)
1184 for (r = 0; r <= R -
p; ++r)
1202 "Method only implemented if BasisType is identical"
1203 "in x and y directions");
1206 "Method only implemented for Modified_A BasisType"
1207 "(x and y direction) and ModifiedPyr_C BasisType (z "
1210 int p,
q, r,
P = 0, Q = 0, idx = 0;
1212 int order0 =
m_base[0]->GetNumModes();
1213 int order1 =
m_base[1]->GetNumModes();
1214 int order2 =
m_base[2]->GetNumModes();
1233 ASSERTL0(
false,
"fid must be between 0 and 4");
1236 if (maparray.size() !=
P * Q)
1247 for (
q = 0;
q < Q; ++
q)
1249 for (
p = 0;
p <
P; ++
p)
1257 for (
p = 0;
p <
P; ++
p)
1259 for (r = 0; r < Q -
p; ++r)
1261 maparray[idx++] =
GetMode(
p, 0, r);
1267 maparray[idx++] =
GetMode(1, 0, 0);
1268 maparray[idx++] =
GetMode(0, 0, 1);
1269 for (r = 1; r < Q - 1; ++r)
1271 maparray[idx++] =
GetMode(1, 0, r);
1274 for (
q = 1;
q <
P; ++
q)
1276 for (r = 0; r < Q -
q; ++r)
1278 maparray[idx++] =
GetMode(1,
q, r);
1284 maparray[idx++] =
GetMode(0, 1, 0);
1285 maparray[idx++] =
GetMode(0, 0, 1);
1286 for (r = 1; r < Q - 1; ++r)
1288 maparray[idx++] =
GetMode(0, 1, r);
1291 for (
p = 1;
p <
P; ++
p)
1293 for (r = 0; r < Q -
p; ++r)
1295 maparray[idx++] =
GetMode(
p, 1, r);
1301 for (
q = 0;
q <
P; ++
q)
1303 for (r = 0; r < Q -
q; ++r)
1305 maparray[idx++] =
GetMode(0,
q, r);
1311 ASSERTL0(
false,
"Face to element map unavailable.");
1321 "Method only implemented if BasisType is identical"
1322 "in x and y directions");
1325 "Method only implemented for Modified_A BasisType"
1326 "(x and y direction) and ModifiedPyr_C BasisType (z "
1329 int i, j, k,
p, r, nFaceCoeffs;
1330 int nummodesA = 0, nummodesB = 0;
1332 int order0 =
m_base[0]->GetNumModes();
1333 int order1 =
m_base[1]->GetNumModes();
1334 int order2 =
m_base[2]->GetNumModes();
1353 ASSERTL0(
false,
"fid must be between 0 and 4");
1364 nFaceCoeffs =
P * (2 * Q -
P + 1) / 2;
1368 nFaceCoeffs =
P * Q;
1372 if (maparray.size() != nFaceCoeffs)
1377 if (signarray.size() != nFaceCoeffs)
1383 fill(signarray.get(), signarray.get() + nFaceCoeffs, 1);
1393 int minPA = min(nummodesA,
P);
1394 int minQB = min(nummodesB, Q);
1396 for (j = 0; j < minPA; ++j)
1399 for (k = 0; k < minQB - j; ++k, ++cnt)
1401 maparray[idx++] = cnt;
1404 cnt += nummodesB - minQB;
1406 for (k = nummodesB - j; k < Q - j; ++k)
1408 signarray[idx] = 0.0;
1409 maparray[idx++] = maparray[0];
1414 for (j = minPA; j < nummodesA; ++j)
1417 for (k = 0; k < minQB-j; ++k, ++cnt)
1419 maparray[idx++] = cnt;
1422 cnt += nummodesB-minQB;
1424 for (k = nummodesB-j; k < Q-j; ++k)
1426 signarray[idx] = 0.0;
1427 maparray[idx++] = maparray[0];
1431 for (j = nummodesA; j <
P; ++j)
1433 for (k = 0; k < Q - j; ++k)
1435 signarray[idx] = 0.0;
1436 maparray[idx++] = maparray[0];
1444 swap(maparray[0], maparray[Q]);
1445 for (i = 1; i < Q - 1; ++i)
1447 swap(maparray[i + 1], maparray[Q + i]);
1451 for (
p = 0;
p <
P; ++
p)
1453 for (r = 0; r < Q -
p; ++r, idx++)
1457 signarray[idx] =
p % 2 ? -1 : 1;
1470 for (i = 0; i < Q; i++)
1472 for (j = 0; j <
P; j++)
1476 arrayindx[i *
P + j] = i *
P + j;
1480 arrayindx[i *
P + j] = j * Q + i;
1487 for (j = 0; j <
P; ++j)
1490 for (k = 0; k < Q; k++)
1492 maparray[arrayindx[j + k *
P]] = j + k * nummodesA;
1495 for (k = nummodesB; k < Q; ++k)
1497 signarray[arrayindx[j + k *
P]] = 0.0;
1498 maparray[arrayindx[j + k *
P]] = maparray[0];
1502 for (j = nummodesA; j <
P; ++j)
1504 for (k = 0; k < Q; ++k)
1506 signarray[arrayindx[j + k *
P]] = 0.0;
1507 maparray[arrayindx[j + k *
P]] = maparray[0];
1522 for (i = 3; i < Q; i += 2)
1524 for (j = 0; j <
P; j++)
1526 signarray[arrayindx[i *
P + j]] *= -1;
1530 for (i = 0; i <
P; i++)
1532 swap(maparray[i], maparray[i +
P]);
1533 swap(signarray[i], signarray[i +
P]);
1538 for (i = 0; i < Q; i++)
1540 for (j = 3; j <
P; j += 2)
1542 signarray[arrayindx[i *
P + j]] *= -1;
1546 for (i = 0; i < Q; i++)
1548 swap(maparray[i], maparray[i + Q]);
1549 swap(signarray[i], signarray[i + Q]);
1561 for (i = 0; i < Q; i++)
1563 for (j = 3; j <
P; j += 2)
1565 signarray[arrayindx[i *
P + j]] *= -1;
1569 for (i = 0; i < Q; i++)
1571 swap(maparray[i *
P], maparray[i *
P + 1]);
1572 swap(signarray[i *
P], signarray[i *
P + 1]);
1577 for (i = 3; i < Q; i += 2)
1579 for (j = 0; j <
P; j++)
1581 signarray[arrayindx[i *
P + j]] *= -1;
1585 for (i = 0; i <
P; i++)
1587 swap(maparray[i * Q], maparray[i * Q + 1]);
1588 swap(signarray[i * Q], signarray[i * Q + 1]);
1601 const int P =
m_base[0]->GetNumModes() - 1;
1602 const int Q =
m_base[1]->GetNumModes() - 1;
1603 const int R =
m_base[2]->GetNumModes() - 1;
1606 if (maparray.size() != nEdgeIntCoeffs)
1611 if (signarray.size() != nEdgeIntCoeffs)
1617 fill(signarray.get(), signarray.get() + nEdgeIntCoeffs, 1);
1627 for (i = 2; i <=
P; ++i)
1629 maparray[i - 2] =
GetMode(i, 0, 0);
1634 for (i = 2; i <= Q; ++i)
1636 maparray[i - 2] =
GetMode(1, i, 0);
1640 for (i = 2; i <=
P; ++i)
1642 maparray[i - 2] =
GetMode(i, 1, 0);
1647 for (i = 2; i <= Q; ++i)
1649 maparray[i - 2] =
GetMode(0, i, 0);
1653 for (i = 2; i <= R; ++i)
1655 maparray[i - 2] =
GetMode(0, 0, i);
1660 for (i = 1; i <= R - 1; ++i)
1662 maparray[i - 1] =
GetMode(1, 0, i);
1666 for (i = 1; i <= R - 1; ++i)
1668 maparray[i - 1] =
GetMode(1, 1, i);
1673 for (i = 1; i <= R - 1; ++i)
1675 maparray[i - 1] =
GetMode(0, 1, i);
1679 ASSERTL0(
false,
"Edge not defined.");
1685 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1696 const int P =
m_base[0]->GetNumModes() - 1;
1697 const int Q =
m_base[1]->GetNumModes() - 1;
1698 const int R =
m_base[2]->GetNumModes() - 1;
1700 int p,
q, r, idx = 0;
1705 if (maparray.size() != nFaceIntCoeffs)
1710 if (signarray.size() != nFaceIntCoeffs)
1716 fill(signarray.get(), signarray.get() + nFaceIntCoeffs, 1);
1727 for (i = 0; i < nummodesB; i++)
1729 for (j = 0; j < nummodesA; j++)
1733 arrayindx[i * nummodesA + j] = i * nummodesA + j;
1737 arrayindx[i * nummodesA + j] = j * nummodesB + i;
1746 for (
q = 2;
q <= Q; ++
q)
1748 for (
p = 2;
p <=
P; ++
p)
1750 maparray[arrayindx[(
q - 2) * nummodesA + (
p - 2)]] =
1756 for (
p = 2;
p <=
P; ++
p)
1758 for (r = 1; r <= R -
p; ++r)
1760 if ((
int)faceOrient == 7)
1762 signarray[idx] =
p % 2 ? -1 : 1;
1764 maparray[idx++] =
GetMode(
p, 0, r);
1769 for (
q = 2;
q <= Q; ++
q)
1771 for (r = 1; r <= R -
q; ++r)
1773 if ((
int)faceOrient == 7)
1775 signarray[idx] =
q % 2 ? -1 : 1;
1777 maparray[idx++] =
GetMode(1,
q, r);
1783 for (
p = 2;
p <=
P; ++
p)
1785 for (r = 1; r <= R -
p; ++r)
1787 if ((
int)faceOrient == 7)
1789 signarray[idx] =
p % 2 ? -1 : 1;
1791 maparray[idx++] =
GetMode(
p, 1, r);
1797 for (
q = 2;
q <= Q; ++
q)
1799 for (r = 1; r <= R -
q; ++r)
1801 if ((
int)faceOrient == 7)
1803 signarray[idx] =
q % 2 ? -1 : 1;
1805 maparray[idx++] =
GetMode(0,
q, r);
1810 ASSERTL0(
false,
"Face interior map not available.");
1820 if (faceOrient == 6 || faceOrient == 8 || faceOrient == 11 ||
1825 for (i = 1; i < nummodesB; i += 2)
1827 for (j = 0; j < nummodesA; j++)
1829 signarray[arrayindx[i * nummodesA + j]] *= -1;
1835 for (i = 0; i < nummodesB; i++)
1837 for (j = 1; j < nummodesA; j += 2)
1839 signarray[arrayindx[i * nummodesA + j]] *= -1;
1845 if (faceOrient == 7 || faceOrient == 8 || faceOrient == 10 ||
1850 for (i = 0; i < nummodesB; i++)
1852 for (j = 1; j < nummodesA; j += 2)
1854 signarray[arrayindx[i * nummodesA + j]] *= -1;
1860 for (i = 1; i < nummodesB; i += 2)
1862 for (j = 0; j < nummodesA; j++)
1864 signarray[arrayindx[i * nummodesA + j]] *= -1;
1911 const int Q =
m_base[1]->GetNumModes() - 1;
1912 const int R =
m_base[2]->GetNumModes() - 1;
1918 for (i = 0; i < I; ++i)
1924 cnt += (R + 1 - i) * (Q + 1) - l * (l + 1) / 2;
1928 l = max(0, J - 1 - I);
1929 cnt += (R + 1 - I) * J - l * (l + 1) / 2;
1943 int nquad0 =
m_base[0]->GetNumPoints();
1944 int nquad1 =
m_base[1]->GetNumPoints();
1945 int nquad2 =
m_base[2]->GetNumPoints();
1954 for (i = 0; i < nquad1 * nquad2; ++i)
1956 Vmath::Vmul(nquad0, inarray.get() + i * nquad0, 1, w0.get(), 1,
1957 outarray.get() + i * nquad0, 1);
1961 for (j = 0; j < nquad2; ++j)
1963 for (i = 0; i < nquad1; ++i)
1966 &outarray[0] + i * nquad0 + j * nquad0 * nquad1, 1);
1976 case LibUtilities::eGaussRadauMAlpha2Beta0:
1977 for (i = 0; i < nquad2; ++i)
1980 &outarray[0] + i * nquad0 * nquad1, 1);
1985 for (i = 0; i < nquad2; ++i)
1988 0.25 * (1 - z2[i]) * (1 - z2[i]) * w2[i],
1989 &outarray[0] + i * nquad0 * nquad1, 1);
1999 int qa =
m_base[0]->GetNumPoints();
2000 int qb =
m_base[1]->GetNumPoints();
2001 int qc =
m_base[2]->GetNumPoints();
2002 int nmodes_a =
m_base[0]->GetNumModes();
2003 int nmodes_b =
m_base[1]->GetNumModes();
2004 int nmodes_c =
m_base[2]->GetNumModes();
2016 int i, j, k, cnt = 0;
2019 OrthoExp.
FwdTrans(array, orthocoeffs);
2029 for (i = 0; i < nmodes_a; ++i)
2031 for (j = 0; j < nmodes_b; ++j)
2033 int maxij = max(i, j);
2035 pow((1.0 * i) / (nmodes_a - 1), cutoff * nmodes_a),
2036 pow((1.0 * j) / (nmodes_b - 1), cutoff * nmodes_b));
2038 for (k = 0; k < nmodes_c - maxij; ++k)
2041 std::max(fac1, pow((1.0 * k) / (nmodes_c - 1),
2042 cutoff * nmodes_c));
2044 orthocoeffs[cnt] *= SvvDiffCoeff * fac;
2060 max_abc = max(max_abc, 0);
2063 for (i = 0; i < nmodes_a; ++i)
2065 for (j = 0; j < nmodes_b; ++j)
2067 int maxij = max(i, j);
2069 for (k = 0; k < nmodes_c - maxij; ++k)
2071 int maxijk = max(maxij, k);
2093 int cutoff_a = (int)(SVVCutOff * nmodes_a);
2094 int cutoff_b = (int)(SVVCutOff * nmodes_b);
2095 int cutoff_c = (int)(SVVCutOff * nmodes_c);
2099 int nmodes = min(min(nmodes_a, nmodes_b), nmodes_c);
2100 NekDouble cutoff = min(min(cutoff_a, cutoff_b), cutoff_c);
2102 for (i = 0; i < nmodes_a; ++i)
2104 for (j = 0; j < nmodes_b; ++j)
2106 int maxij = max(i, j);
2107 for (k = 0; k < nmodes_c - maxij; ++k)
2109 if (j + k >= cutoff || i + k >= cutoff)
2113 exp(-(i + k - nmodes) * (i + k - nmodes) /
2114 ((
NekDouble)((i + k - cutoff + epsilon) *
2115 (i + k - cutoff + epsilon)))) *
2116 exp(-(j - nmodes) * (j - nmodes) /
2118 (j - cutoff + epsilon)))));
2122 orthocoeffs[cnt] *= 0.0;
2131 OrthoExp.
BwdTrans(orthocoeffs, array);
2138 int nquad0 =
m_base[0]->GetNumPoints();
2139 int nquad1 =
m_base[1]->GetNumPoints();
2140 int nquad2 =
m_base[2]->GetNumPoints();
2141 int nqtot = nquad0 * nquad1 * nquad2;
2142 int nmodes0 =
m_base[0]->GetNumModes();
2143 int nmodes1 =
m_base[1]->GetNumModes();
2144 int nmodes2 =
m_base[2]->GetNumModes();
2145 int numMax = nmodes0;
2166 bortho0, bortho1, bortho2);
2169 OrthoPyrExp->FwdTrans(phys_tmp, coeff);
2172 for (u = 0; u < numMin; ++u)
2174 for (i = 0; i < numMin; ++i)
2177 int maxui = max(u, i);
2179 tmp2 = coeff_tmp1 + cnt, 1);
2180 cnt += nmodes2 - maxui;
2183 for (i = numMin; i < nmodes1; ++i)
2185 int maxui = max(u, i);
2186 cnt += numMax - maxui;
2190 OrthoPyrExp->BwdTrans(coeff_tmp1, phys_tmp);
#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.
BasisType GetBasisType() const
Return type of expansion basis.
int GetNumModes() const
Returns the order of the basis.
Defines a specification for a set of points.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
NekDouble BaryTensorDeriv(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d1, Array< OneD, NekDouble > &outarray_d2, Array< OneD, NekDouble > &outarray_d3)
Calculate the 3D derivative in the local tensor/collapsed coordinate at the physical points.
The base class for all shapes.
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
int NumBndryCoeffs(void) const
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
const Array< OneD, const LibUtilities::BasisSharedPtr > & GetBase() const
This function gets the shared point to basis.
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
NekDouble PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
This function evaluates the expansion at a single (arbitrary) point of the domain.
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Backward transformation from coefficient space to physical space.
int GetTraceNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th trace.
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Forward transformation from physical space to coefficient space.
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
void PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Array< OneD, LibUtilities::BasisSharedPtr > m_base
NekDouble GetConstFactor(const ConstFactorType &factor) const
bool ConstFactorExists(const ConstFactorType &factor) const
void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override
int v_NumBndryCoeffs() const override
NekDouble v_PhysEvaluateBasis(const Array< OneD, const NekDouble > &coords, int mode) final
void v_ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
int v_NumDGBndryCoeffs() const override
void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Forward transform from physical quadrature space stored in inarray and evaluate the expansion coeffic...
void v_GetElmtTraceToTraceMap(const unsigned int fid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation faceOrient, int P, int Q) override
int v_GetTraceNcoeffs(const int i) const override
void v_LocCollapsedToLocCoord(const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi) override
void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey) override
int v_GetEdgeNcoeffs(const int i) const override
void v_BwdTrans_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2) override
int v_GetNtraces() const override
void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray) override
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
int v_GetNedges() const override
int v_CalcNumberOfCoefficients(const std::vector< unsigned int > &nummodes, int &modes_offset) override
int GetMode(int I, int J, int K)
Compute the mode number in the expansion for a particular tensorial combination.
int v_GetVertexMap(int localVertexId, bool useCoeffPacking=false) override
int v_GetTraceNumPoints(const int i) const override
void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Backward transformation is evaluated at the quadrature points.
void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray) override
void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Inner product of inarray over region with respect to the expansion basis m_base[0]->GetBdata(),...
void v_GetInteriorMap(Array< OneD, unsigned int > &outarray) override
void v_GetEdgeInteriorToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2) override
void v_LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta) override
DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey) override
int v_GetNverts() const override
void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_GetTraceInteriorToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2) override
void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2) override
void v_StdPhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2) override
const LibUtilities::BasisKey v_GetTraceBasisKey(const int i, const int k) const override
LibUtilities::ShapeType v_DetShapeType() const override
void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2) override
Calculate the derivative of the physical points.
void v_GetTraceNumModes(const int fid, int &numModes0, int &numModes1, Orientation faceOrient=eDir1FwdDir1_Dir2FwdDir2) override
NekDouble v_PhysEvaluate(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
void v_GetTraceCoeffMap(const unsigned int fid, Array< OneD, unsigned int > &maparray) override
void v_GetCoords(Array< OneD, NekDouble > &xi_x, Array< OneD, NekDouble > &xi_y, Array< OneD, NekDouble > &xi_z) override
int v_GetTraceIntNcoeffs(const int i) const override
static void Dgemv(const char &trans, const int &m, const int &n, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = alpha A x plus beta y where A[m x n].
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
static double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = .
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...
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.
int getNumberOfCoefficients(int Na, int Nb, int Nc)
int getNumberOfBndCoefficients(int Na, int Nb, int Nc)
int getNumberOfCoefficients(int Na)
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.
@ P
Monomial polynomials .
@ eOrtho_A
Principle Orthogonal Functions .
@ eGLL_Lagrange
Lagrange for SEM basis .
@ eModifiedPyr_C
Principle Modified Functions.
@ eModified_A
Principle Modified Functions .
@ eOrthoPyr_C
Principle Orthogonal Functions .
static const NekDouble kNekZeroTol
LibUtilities::BasisKey EvaluateTriFaceBasisKey(const int facedir, const LibUtilities::BasisType faceDirBasisType, const int numpoints, const int nummodes)
@ eFactorSVVDGKerDiffCoeff
@ eFactorSVVPowerKerDiffCoeff
const int kSVVDGFiltermodesmin
std::shared_ptr< StdPyrExp > StdPyrExpSharedPtr
const int kSVVDGFiltermodesmax
LibUtilities::BasisKey EvaluateQuadFaceBasisKey(const int facedir, const LibUtilities::BasisType faceDirBasisType, const int numpoints, const int nummodes)
const NekDouble kSVVDGFilter[9][11]
@ eDir1BwdDir2_Dir2BwdDir1
@ eDir1BwdDir1_Dir2BwdDir2
@ eDir1BwdDir2_Dir2FwdDir1
@ eDir1FwdDir1_Dir2BwdDir2
@ eDir1BwdDir1_Dir2FwdDir2
@ eDir1FwdDir2_Dir2FwdDir1
@ eDir1FwdDir2_Dir2BwdDir1
std::vector< double > q(NPUPPER *NPUPPER)
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.
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.
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)