46 :
StdExpansion(LibUtilities::StdPrismData::getNumberOfCoefficients(
47 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
49 StdExpansion3D(LibUtilities::StdPrismData::getNumberOfCoefficients(
50 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
54 "order in 'a' direction is higher than order in 'c' direction");
80 int Qx =
m_base[0]->GetNumPoints();
81 int Qy =
m_base[1]->GetNumPoints();
82 int Qz =
m_base[2]->GetNumPoints();
83 int Qtot = Qx * Qy * Qz;
93 bool Do_1 = (out_dxi1.size() > 0) ?
true :
false;
94 bool Do_3 = (out_dxi3.size() > 0) ?
true :
false;
113 for (k = 0; k < Qz; ++k)
116 &dEta_bar1[0] + k * Qx * Qy, 1,
117 &out_dxi1[0] + k * Qx * Qy, 1);
124 for (k = 0; k < Qz; ++k)
127 &dEta_bar1[0] + k * Qx * Qy, 1,
128 &dEta_bar1[0] + k * Qx * Qy, 1);
132 for (i = 0; i < Qx; ++i)
134 Vmath::Svtvp(Qz * Qy, 1.0 + eta_x[i], &dEta_bar1[0] + i, Qx,
135 &out_dxi3[0] + i, Qx, &out_dxi3[0] + i, Qx);
169 ASSERTL1(
false,
"input dir is out of range");
218 "Basis[1] is not a general tensor type");
222 "Basis[2] is not a general tensor type");
224 if (
m_base[0]->Collocation() &&
m_base[1]->Collocation() &&
229 inarray, 1, outarray, 1);
240 int nquad1 =
m_base[1]->GetNumPoints();
241 int nquad2 =
m_base[2]->GetNumPoints();
242 int order0 =
m_base[0]->GetNumModes();
243 int order1 =
m_base[1]->GetNumModes();
246 nquad1 * nquad2 * order0);
249 m_base[2]->GetBdata(), inarray, outarray, wsp,
true,
259 [[maybe_unused]]
bool doCheckCollDir0,
260 [[maybe_unused]]
bool doCheckCollDir1,
261 [[maybe_unused]]
bool doCheckCollDir2)
264 int nquad0 =
m_base[0]->GetNumPoints();
265 int nquad1 =
m_base[1]->GetNumPoints();
266 int nquad2 =
m_base[2]->GetNumPoints();
267 int nummodes0 =
m_base[0]->GetNumModes();
268 int nummodes1 =
m_base[1]->GetNumModes();
269 int nummodes2 =
m_base[2]->GetNumModes();
273 for (i = mode = 0; i < nummodes0; ++i)
275 Blas::Dgemm(
'N',
'N', nquad2, nummodes1, nummodes2 - i, 1.0,
276 base2.data() + mode * nquad2, nquad2,
277 inarray.data() + mode * nummodes1, nummodes2 - i, 0.0,
278 tmp0.data() + i * nquad2 * nummodes1, nquad2);
279 mode += nummodes2 - i;
284 for (i = 0; i < nummodes1; i++)
287 base2.data() + nquad2, 1,
288 tmp0.data() + nquad2 * (nummodes1 + i), 1);
292 for (i = 0; i < nummodes0; i++)
294 Blas::Dgemm(
'N',
'T', nquad1, nquad2, nummodes1, 1.0, base1.data(),
295 nquad1, tmp0.data() + i * nquad2 * nummodes1, nquad2, 0.0,
296 tmp1.data() + i * nquad2 * nquad1, nquad1);
299 Blas::Dgemm(
'N',
'T', nquad0, nquad2 * nquad1, nummodes0, 1.0, base0.data(),
300 nquad0, tmp1.data(), nquad2 * nquad1, 0.0, outarray.data(),
328 out = (*matsys) * in;
366 "Basis[1] is not a general tensor type");
370 "Basis[2] is not a general tensor type");
372 if (
m_base[0]->Collocation() &&
m_base[1]->Collocation())
386 int nquad1 =
m_base[1]->GetNumPoints();
387 int nquad2 =
m_base[2]->GetNumPoints();
388 int order0 =
m_base[0]->GetNumModes();
389 int order1 =
m_base[1]->GetNumModes();
393 if (multiplybyweights)
400 tmp, outarray, wsp,
true,
true,
true);
406 inarray, outarray, wsp,
true,
true,
true);
416 [[maybe_unused]]
bool doCheckCollDir0,
417 [[maybe_unused]]
bool doCheckCollDir1,
418 [[maybe_unused]]
bool doCheckCollDir2)
422 const int nquad0 =
m_base[0]->GetNumPoints();
423 const int nquad1 =
m_base[1]->GetNumPoints();
424 const int nquad2 =
m_base[2]->GetNumPoints();
425 const int order0 =
m_base[0]->GetNumModes();
426 const int order1 =
m_base[1]->GetNumModes();
427 const int order2 =
m_base[2]->GetNumModes();
431 ASSERTL1(wsp.size() >= nquad1 * nquad2 * order0 + nquad2 * order0 * order1,
432 "Insufficient workspace size");
438 Blas::Dgemm(
'T',
'N', nquad1 * nquad2, order0, nquad0, 1.0, inarray.data(),
439 nquad0, base0.data(), nquad0, 0.0, tmp0.data(),
443 Blas::Dgemm(
'T',
'N', nquad2 * order0, order1, nquad1, 1.0, tmp0.data(),
444 nquad1, base1.data(), nquad1, 0.0, tmp1.data(),
448 for (mode = i = 0; i < order0; ++i)
450 Blas::Dgemm(
'T',
'N', order2 - i, order1, nquad2, 1.0,
451 base2.data() + mode * nquad2, nquad2,
452 tmp1.data() + i * nquad2, nquad2 * order0, 0.0,
453 outarray.data() + mode * order1, order2 - i);
461 for (i = 0; i < order1; ++i)
466 tmp1.data() + i * order0 * nquad2 + nquad2, 1);
486 ASSERTL0(dir >= 0 && dir <= 2,
"input dir is out of range");
489 int order0 =
m_base[0]->GetNumModes();
490 int order1 =
m_base[1]->GetNumModes();
491 int nquad0 =
m_base[0]->GetNumPoints();
492 int nquad1 =
m_base[1]->GetNumPoints();
493 int nquad2 =
m_base[2]->GetNumPoints();
503 for (i = 0; i < nquad0; ++i)
505 gfac0[i] = 0.5 * (1 + z0[i]);
509 for (i = 0; i < nquad2; ++i)
511 gfac2[i] = 2.0 / (1 - z2[i]);
517 for (i = 0; i < nquad2; ++i)
520 &inarray[0] + i * nquad0 * nquad1, 1,
521 &tmp0[0] + i * nquad0 * nquad1, 1);
532 m_base[2]->GetBdata(), tmp0, outarray, wsp,
true,
true,
true);
540 m_base[2]->GetBdata(), tmp0, outarray, wsp,
true,
true,
true);
549 for (i = 0; i < nquad1 * nquad2; ++i)
551 Vmath::Vmul(nquad0, &gfac0[0], 1, &tmp0[0] + i * nquad0, 1,
552 &tmp0[0] + i * nquad0, 1);
557 m_base[2]->GetBdata(), tmp0, tmp1, wsp,
true,
true,
true);
562 m_base[2]->GetDbdata(), tmp0, outarray, wsp,
true,
true,
true);
592 eta[0] = 2.0 * (1.0 + xi[0]) / d2 - 1.0;
598 xi[0] = (1.0 + eta[0]) * (1.0 - eta[2]) * 0.5 - 1.0;
615 for (
int k = 0; k < Qz; ++k)
617 for (
int j = 0; j < Qy; ++j)
619 for (
int i = 0; i < Qx; ++i)
621 int s = i + Qx * (j + Qy * k);
622 xi_x[s] = (1.0 - eta_z[k]) * (1.0 + etaBar_x[i]) / 2.0 - 1.0;
636 const int nm1 =
m_base[1]->GetNumModes();
637 const int nm2 =
m_base[2]->GetNumModes();
638 const int b = 2 * nm2 + 1;
640 const int mode0 = floor(0.5 * (b -
sqrt(b * b - 8.0 * mode / nm1)));
642 mode - nm1 * (mode0 * (nm2 - 1) + 1 - (mode0 - 2) * (mode0 - 1) / 2);
643 const int mode1 = tmp / (nm2 - mode0);
644 const int mode2 = tmp % (nm2 - mode0);
646 if (mode0 == 0 && mode2 == 1 &&
650 return StdExpansion::BaryEvaluateBasis<1>(coll[1], mode1) *
651 StdExpansion::BaryEvaluateBasis<2>(coll[2], mode2);
655 return StdExpansion::BaryEvaluateBasis<0>(coll[0], mode0) *
656 StdExpansion::BaryEvaluateBasis<1>(coll[1], mode1) *
657 StdExpansion::BaryEvaluateBasis<2>(coll[2], mode2);
664 std::array<NekDouble, 3> &firstOrderDerivs)
673 if ((1 - coll[2]) < 1e-5)
677 EphysDeriv2(totPoints);
678 PhysDeriv(inarray, EphysDeriv0, EphysDeriv1, EphysDeriv2);
681 I[0] =
GetBase()[0]->GetI(coll);
682 I[1] =
GetBase()[1]->GetI(coll + 1);
683 I[2] =
GetBase()[2]->GetI(coll + 2);
693 NekDouble dEta_bar1 = firstOrderDerivs[0];
696 firstOrderDerivs[0] = fac * dEta_bar1;
699 fac = 1.0 / (1.0 - coll[2]);
700 dEta_bar1 = fac * dEta_bar1;
704 firstOrderDerivs[2] += fac * dEta_bar1;
719 int nummodes[3] = {
m_base[0]->GetNumModes(),
m_base[1]->GetNumModes(),
720 m_base[2]->GetNumModes()};
726 numModes0 = nummodes[0];
727 numModes1 = nummodes[1];
734 numModes0 = nummodes[1];
735 numModes1 = nummodes[2];
742 numModes0 = nummodes[0];
743 numModes1 = nummodes[2];
750 std::swap(numModes0, numModes1);
756 ASSERTL2(i >= 0 && i <= 8,
"edge id is out of range");
758 if (i == 0 || i == 2)
762 else if (i == 1 || i == 3 || i == 8)
804 "BasisType is not a boundary interior form");
807 "BasisType is not a boundary interior form");
810 "BasisType is not a boundary interior form");
812 int P =
m_base[0]->GetNumModes();
813 int Q =
m_base[1]->GetNumModes();
814 int R =
m_base[2]->GetNumModes();
823 "BasisType is not a boundary interior form");
826 "BasisType is not a boundary interior form");
829 "BasisType is not a boundary interior form");
831 int P =
m_base[0]->GetNumModes() - 1;
832 int Q =
m_base[1]->GetNumModes() - 1;
833 int R =
m_base[2]->GetNumModes() - 1;
835 return (
P + 1) * (Q + 1)
836 + 2 * (Q + 1) * (R + 1)
837 + 2 * (R + 1) +
P * (1 + 2 * R -
P);
842 ASSERTL2(i >= 0 && i <= 4,
"face id is out of range");
847 else if (i == 1 || i == 3)
850 return Q + 1 + (
P * (1 + 2 * Q -
P)) / 2;
860 ASSERTL2(i >= 0 && i <= 4,
"face id is out of range");
870 else if (i == 1 || i == 3)
872 return Pi * (2 * Ri - Pi - 1) / 2;
882 ASSERTL2(i >= 0 && i <= 4,
"face id is out of range");
886 return m_base[0]->GetNumPoints() *
m_base[1]->GetNumPoints();
888 else if (i == 1 || i == 3)
890 return m_base[0]->GetNumPoints() *
m_base[2]->GetNumPoints();
894 return m_base[1]->GetNumPoints() *
m_base[2]->GetNumPoints();
901 ASSERTL2(i >= 0 && i <= 4,
"face id is out of range");
902 ASSERTL2(j == 0 || j == 1,
"face direction is out of range");
906 return m_base[j]->GetPointsKey();
908 else if (i == 1 || i == 3)
910 return m_base[2 * j]->GetPointsKey();
914 return m_base[j + 1]->GetPointsKey();
922 ASSERTL2(i >= 0 && i <= 4,
"face id is out of range");
923 ASSERTL2(k >= 0 && k <= 1,
"basis key id is out of range");
949 const std::vector<unsigned int> &nummodes,
int &modes_offset)
952 nummodes[modes_offset], nummodes[modes_offset + 1],
953 nummodes[modes_offset + 2]);
975 "Mapping not defined for this type of basis");
979 if (useCoeffPacking ==
true)
1002 ASSERTL0(
false,
"local vertex id must be between 0 and 5");
1028 ASSERTL0(
false,
"local vertex id must be between 0 and 5");
1039 "BasisType is not a boundary interior form");
1042 "BasisType is not a boundary interior form");
1045 "BasisType is not a boundary interior form");
1047 int P =
m_base[0]->GetNumModes() - 1, p;
1048 int Q =
m_base[1]->GetNumModes() - 1, q;
1049 int R =
m_base[2]->GetNumModes() - 1, r;
1053 if (outarray.size() != nIntCoeffs)
1061 for (p = 2; p <=
P; ++p)
1063 for (q = 2; q <= Q; ++q)
1065 for (r = 1; r <= R - p; ++r)
1067 outarray[idx++] =
GetMode(p, q, r);
1077 "BasisType is not a boundary interior form");
1080 "BasisType is not a boundary interior form");
1083 "BasisType is not a boundary interior form");
1085 int P =
m_base[0]->GetNumModes() - 1, p;
1086 int Q =
m_base[1]->GetNumModes() - 1, q;
1087 int R =
m_base[2]->GetNumModes() - 1, r;
1092 if (maparray.size() != nBnd)
1098 for (p = 0; p <=
P; ++p)
1103 for (q = 0; q <= Q; ++q)
1105 for (r = 0; r <= R - p; ++r)
1107 maparray[idx++] =
GetMode(p, q, r);
1115 for (q = 0; q <= Q; ++q)
1119 for (r = 0; r <= R - p; ++r)
1121 maparray[idx++] =
GetMode(p, q, r);
1126 maparray[idx++] =
GetMode(p, q, 0);
1137 "Method only implemented if BasisType is identical"
1138 "in x and y directions");
1141 "Method only implemented for Modified_A BasisType"
1142 "(x and y direction) and Modified_B BasisType (z "
1144 int p, q, r, idx = 0;
1151 Q =
m_base[1]->GetNumModes();
1156 Q =
m_base[2]->GetNumModes();
1161 Q =
m_base[2]->GetNumModes();
1164 ASSERTL0(
false,
"fid must be between 0 and 4");
1167 if (maparray.size() !=
P * Q)
1177 for (q = 0; q < Q; ++q)
1179 for (p = 0; p <
P; ++p)
1181 maparray[q *
P + p] =
GetMode(p, q, 0);
1186 for (p = 0; p <
P; ++p)
1188 for (r = 0; r < Q - p; ++r)
1190 maparray[idx++] =
GetMode(p, 0, r);
1195 for (q = 0; q <
P; ++q)
1197 maparray[q] =
GetMode(1, q, 0);
1199 for (q = 0; q <
P; ++q)
1201 maparray[
P + q] =
GetMode(0, q, 1);
1203 for (r = 1; r < Q - 1; ++r)
1205 for (q = 0; q <
P; ++q)
1207 maparray[(r + 1) *
P + q] =
GetMode(1, q, r);
1212 for (p = 0; p <
P; ++p)
1214 for (r = 0; r < Q - p; ++r)
1216 maparray[idx++] =
GetMode(p, 1, r);
1221 for (r = 0; r < Q; ++r)
1223 for (q = 0; q <
P; ++q)
1225 maparray[r *
P + q] =
GetMode(0, q, r);
1230 ASSERTL0(
false,
"Face to element map unavailable.");
1240 "Method only implemented if BasisType is identical"
1241 "in x and y directions");
1244 "Method only implemented for Modified_A BasisType"
1245 "(x and y direction) and Modified_B BasisType (z "
1248 int i, j, k, p, r, nFaceCoeffs, idx = 0;
1249 int nummodesA = 0, nummodesB = 0;
1254 nummodesA =
m_base[0]->GetNumModes();
1255 nummodesB =
m_base[1]->GetNumModes();
1259 nummodesA =
m_base[0]->GetNumModes();
1260 nummodesB =
m_base[2]->GetNumModes();
1264 nummodesA =
m_base[1]->GetNumModes();
1265 nummodesB =
m_base[2]->GetNumModes();
1268 ASSERTL0(
false,
"fid must be between 0 and 4");
1277 else if (fid == 1 || fid == 3)
1279 nFaceCoeffs =
P * (2 * Q -
P + 1) / 2;
1283 nFaceCoeffs =
P * Q;
1287 if (maparray.size() != nFaceCoeffs)
1292 if (signarray.size() != nFaceCoeffs)
1298 fill(signarray.data(), signarray.data() + nFaceCoeffs, 1);
1301 int minPA =
min(nummodesA,
P);
1302 int minQB =
min(nummodesB, Q);
1304 if (fid == 1 || fid == 3)
1311 for (j = 0; j < minPA; ++j)
1314 for (k = 0; k < minQB - j; ++k, ++cnt)
1316 maparray[idx++] = cnt;
1319 cnt += nummodesB - minQB;
1322 for (k = nummodesB - j; k < Q - j; ++k)
1324 signarray[idx] = 0.0;
1325 maparray[idx++] = maparray[0];
1329 for (j = minPA; j < nummodesA; ++j)
1332 for (k = 0; k < minQB-j; ++k, ++cnt)
1334 maparray[idx++] = cnt;
1337 cnt += nummodesB-minQB;
1340 for (k = nummodesB-j; k < Q-j; ++k)
1342 signarray[idx] = 0.0;
1343 maparray[idx++] = maparray[0];
1347 for (j = nummodesA; j <
P; ++j)
1349 for (k = 0; k < Q - j; ++k)
1351 signarray[idx] = 0.0;
1352 maparray[idx++] = maparray[0];
1361 for (p = 0; p <
P; ++p)
1363 for (r = 0; r < Q - p; ++r, idx++)
1367 signarray[idx] = p % 2 ? -1 : 1;
1372 swap(maparray[0], maparray[Q]);
1373 for (i = 1; i < Q - 1; ++i)
1375 swap(maparray[i + 1], maparray[Q + i]);
1385 for (i = 0; i < Q; i++)
1387 for (j = 0; j <
P; j++)
1391 arrayindx[i *
P + j] = i *
P + j;
1395 arrayindx[i *
P + j] = j * Q + i;
1402 for (j = 0; j <
P; ++j)
1405 for (k = 0; k < Q; k++)
1407 maparray[arrayindx[j + k *
P]] = j + k * nummodesA;
1410 for (k = nummodesB; k < Q; ++k)
1412 signarray[arrayindx[j + k *
P]] = 0.0;
1413 maparray[arrayindx[j + k *
P]] = maparray[0];
1417 for (j = nummodesA; j <
P; ++j)
1419 for (k = 0; k < Q; ++k)
1421 signarray[arrayindx[j + k *
P]] = 0.0;
1422 maparray[arrayindx[j + k *
P]] = maparray[0];
1437 for (i = 3; i < Q; i += 2)
1439 for (j = 0; j <
P; j++)
1441 signarray[arrayindx[i *
P + j]] *= -1;
1445 for (i = 0; i <
P; i++)
1447 swap(maparray[i], maparray[i +
P]);
1448 swap(signarray[i], signarray[i +
P]);
1453 for (i = 0; i < Q; i++)
1455 for (j = 3; j <
P; j += 2)
1457 signarray[arrayindx[i *
P + j]] *= -1;
1461 for (i = 0; i < Q; i++)
1463 swap(maparray[i], maparray[i + Q]);
1464 swap(signarray[i], signarray[i + Q]);
1476 for (i = 0; i < Q; i++)
1478 for (j = 3; j <
P; j += 2)
1480 signarray[arrayindx[i *
P + j]] *= -1;
1484 for (i = 0; i < Q; i++)
1486 swap(maparray[i *
P], maparray[i *
P + 1]);
1487 swap(signarray[i *
P], signarray[i *
P + 1]);
1492 for (i = 3; i < Q; i += 2)
1494 for (j = 0; j <
P; j++)
1496 signarray[arrayindx[i *
P + j]] *= -1;
1500 for (i = 0; i <
P; i++)
1502 swap(maparray[i * Q], maparray[i * Q + 1]);
1503 swap(signarray[i * Q], signarray[i * Q + 1]);
1516 const int P =
m_base[0]->GetNumModes() - 1;
1517 const int Q =
m_base[1]->GetNumModes() - 1;
1518 const int R =
m_base[2]->GetNumModes() - 1;
1521 if (maparray.size() != nEdgeIntCoeffs)
1526 if (signarray.size() != nEdgeIntCoeffs)
1532 fill(signarray.data(), signarray.data() + nEdgeIntCoeffs, 1);
1542 for (i = 2; i <=
P; ++i)
1544 maparray[i - 2] =
GetMode(i, 0, 0);
1549 for (i = 2; i <= Q; ++i)
1551 maparray[i - 2] =
GetMode(1, i, 0);
1558 for (i = 2; i <=
P; ++i)
1560 maparray[i - 2] =
GetMode(i, 1, 0);
1567 for (i = 2; i <= Q; ++i)
1569 maparray[i - 2] =
GetMode(0, i, 0);
1574 for (i = 2; i <= R; ++i)
1576 maparray[i - 2] =
GetMode(0, 0, i);
1581 for (i = 1; i <= R - 1; ++i)
1583 maparray[i - 1] =
GetMode(1, 0, i);
1588 for (i = 1; i <= R - 1; ++i)
1590 maparray[i - 1] =
GetMode(1, 1, i);
1595 for (i = 2; i <= R; ++i)
1597 maparray[i - 2] =
GetMode(0, 1, i);
1602 for (i = 2; i <= Q; ++i)
1604 maparray[i - 2] =
GetMode(0, i, 1);
1609 ASSERTL0(
false,
"Edge not defined.");
1615 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1626 const int P =
m_base[0]->GetNumModes() - 1;
1627 const int Q =
m_base[1]->GetNumModes() - 1;
1628 const int R =
m_base[2]->GetNumModes() - 1;
1630 int p, q, r, idx = 0;
1636 if (maparray.size() != nFaceIntCoeffs)
1641 if (signarray.size() != nFaceIntCoeffs)
1647 fill(signarray.data(), signarray.data() + nFaceIntCoeffs, 1);
1653 if (fid != 1 && fid != 3)
1666 for (i = 0; i < nummodesB; i++)
1668 for (j = 0; j < nummodesA; j++)
1672 arrayindx[i * nummodesA + j] = i * nummodesA + j;
1676 arrayindx[i * nummodesA + j] = j * nummodesB + i;
1685 for (q = 2; q <= Q; ++q)
1687 for (p = 2; p <=
P; ++p)
1689 maparray[arrayindx[(q - 2) * nummodesA + (p - 2)]] =
1696 for (p = 2; p <=
P; ++p)
1698 for (r = 1; r <= R - p; ++r)
1702 signarray[idx] = p % 2 ? -1 : 1;
1704 maparray[idx++] =
GetMode(p, 0, r);
1710 for (r = 1; r <= R - 1; ++r)
1712 for (q = 2; q <= Q; ++q)
1714 maparray[arrayindx[(r - 1) * nummodesA + (q - 2)]] =
1721 for (p = 2; p <=
P; ++p)
1723 for (r = 1; r <= R - p; ++r)
1727 signarray[idx] = p % 2 ? -1 : 1;
1729 maparray[idx++] =
GetMode(p, 1, r);
1735 for (r = 2; r <= R; ++r)
1737 for (q = 2; q <= Q; ++q)
1739 maparray[arrayindx[(r - 2) * nummodesA + (q - 2)]] =
1746 ASSERTL0(
false,
"Face interior map not available.");
1751 if (fid == 1 || fid == 3)
1763 for (i = 1; i < nummodesB; i += 2)
1765 for (j = 0; j < nummodesA; j++)
1767 signarray[arrayindx[i * nummodesA + j]] *= -1;
1773 for (i = 0; i < nummodesB; i++)
1775 for (j = 1; j < nummodesA; j += 2)
1777 signarray[arrayindx[i * nummodesA + j]] *= -1;
1790 for (i = 0; i < nummodesB; i++)
1792 for (j = 1; j < nummodesA; j += 2)
1794 signarray[arrayindx[i * nummodesA + j]] *= -1;
1800 for (i = 1; i < nummodesB; i += 2)
1802 for (j = 0; j < nummodesA; j++)
1804 signarray[arrayindx[i * nummodesA + j]] *= -1;
1826 int nq0 =
m_base[0]->GetNumPoints();
1827 int nq1 =
m_base[1]->GetNumPoints();
1828 int nq2 =
m_base[2]->GetNumPoints();
1838 nq =
max(nq0,
max(nq1, nq2));
1851 for (
int i = 0; i < nq; ++i)
1853 for (
int j = 0; j < nq; ++j)
1855 for (
int k = 0; k < nq - i; ++k, ++cnt)
1858 coords[cnt][0] = -1.0 + 2 * k / (
NekDouble)(nq - 1);
1859 coords[cnt][1] = -1.0 + 2 * j / (
NekDouble)(nq - 1);
1860 coords[cnt][2] = -1.0 + 2 * i / (
NekDouble)(nq - 1);
1865 for (
int i = 0; i < neq; ++i)
1869 I[0] =
m_base[0]->GetI(coll);
1870 I[1] =
m_base[1]->GetI(coll + 1);
1871 I[2] =
m_base[2]->GetI(coll + 2);
1875 for (
int k = 0; k < nq2; ++k)
1877 for (
int j = 0; j < nq1; ++j)
1880 fac = (I[1]->GetPtr())[j] * (I[2]->GetPtr())[k];
1884 Mat->GetRawPtr() + k * nq0 * nq1 * neq +
1927 int Q =
m_base[1]->GetNumModes() - 1;
1928 int R =
m_base[2]->GetNumModes() - 1;
1932 (Q + 1) * (p * R + 1 -
1933 (p - 2) * (p - 1) / 2);
1941 int nquad0 =
m_base[0]->GetNumPoints();
1942 int nquad1 =
m_base[1]->GetNumPoints();
1943 int nquad2 =
m_base[2]->GetNumPoints();
1952 for (i = 0; i < nquad1 * nquad2; ++i)
1954 Vmath::Vmul(nquad0, inarray.data() + i * nquad0, 1, w0.data(), 1,
1955 outarray.data() + i * nquad0, 1);
1959 for (j = 0; j < nquad2; ++j)
1961 for (i = 0; i < nquad1; ++i)
1964 &outarray[0] + i * nquad0 + j * nquad0 * nquad1, 1);
1974 case LibUtilities::eGaussRadauMAlpha1Beta0:
1975 for (i = 0; i < nquad2; ++i)
1978 &outarray[0] + i * nquad0 * nquad1, 1);
1983 for (i = 0; i < nquad2; ++i)
1985 Blas::Dscal(nquad0 * nquad1, 0.5 * (1 - z2[i]) * w2[i],
1986 &outarray[0] + i * nquad0 * nquad1, 1);
1996 int qa =
m_base[0]->GetNumPoints();
1997 int qb =
m_base[1]->GetNumPoints();
1998 int qc =
m_base[2]->GetNumPoints();
1999 int nmodes_a =
m_base[0]->GetNumModes();
2000 int nmodes_b =
m_base[1]->GetNumModes();
2001 int nmodes_c =
m_base[2]->GetNumModes();
2013 int i, j, k, cnt = 0;
2016 OrthoExp.
FwdTrans(array, orthocoeffs);
2026 for (i = 0; i < nmodes_a; ++i)
2028 for (j = 0; j < nmodes_b; ++j)
2031 pow((1.0 * i) / (nmodes_a - 1), cutoff * nmodes_a),
2032 pow((1.0 * j) / (nmodes_b - 1), cutoff * nmodes_b));
2034 for (k = 0; k < nmodes_c - i; ++k)
2037 std::max(fac1, pow((1.0 * k) / (nmodes_c - 1),
2038 cutoff * nmodes_c));
2040 orthocoeffs[cnt] *= SvvDiffCoeff * fac;
2056 max_abc =
max(max_abc, 0);
2059 for (i = 0; i < nmodes_a; ++i)
2061 for (j = 0; j < nmodes_b; ++j)
2063 int maxij =
max(i, j);
2065 for (k = 0; k < nmodes_c - i; ++k)
2067 int maxijk =
max(maxij, k);
2089 int cutoff_a = (int)(SVVCutOff * nmodes_a);
2090 int cutoff_b = (int)(SVVCutOff * nmodes_b);
2091 int cutoff_c = (int)(SVVCutOff * nmodes_c);
2095 int nmodes =
min(
min(nmodes_a, nmodes_b), nmodes_c);
2099 for (i = 0; i < nmodes_a; ++i)
2101 for (j = 0; j < nmodes_b; ++j)
2103 for (k = 0; k < nmodes_c - i; ++k)
2105 if (j >= cutoff || i + k >= cutoff)
2109 exp(-(i + k - nmodes) * (i + k - nmodes) /
2110 ((
NekDouble)((i + k - cutoff + epsilon) *
2111 (i + k - cutoff + epsilon)))) *
2112 exp(-(j - nmodes) * (j - nmodes) /
2114 (j - cutoff + epsilon)))));
2118 orthocoeffs[cnt] *= 0.0;
2127 OrthoExp.
BwdTrans(orthocoeffs, array);
2134 int nquad0 =
m_base[0]->GetNumPoints();
2135 int nquad1 =
m_base[1]->GetNumPoints();
2136 int nquad2 =
m_base[2]->GetNumPoints();
2137 int nqtot = nquad0 * nquad1 * nquad2;
2138 int nmodes0 =
m_base[0]->GetNumModes();
2139 int nmodes1 =
m_base[1]->GetNumModes();
2140 int nmodes2 =
m_base[2]->GetNumModes();
2141 int numMax = nmodes0;
2162 bortho0, bortho1, bortho2);
2165 OrthoPrismExp->FwdTrans(phys_tmp, coeff);
2168 for (u = 0; u < numMin; ++u)
2170 for (i = 0; i < numMin; ++i)
2173 tmp2 = coeff_tmp1 + cnt, 1);
2177 for (i = numMin; i < numMax; ++i)
2183 OrthoPrismExp->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.
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 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.
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.
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
Class representing a prismatic element in reference space.
void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey) override
void v_GetElmtTraceToTraceMap(const unsigned int fid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation faceOrient, int P, int Q) override
void v_ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_GetEdgeInteriorToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2) override
LibUtilities::PointsKey v_GetTracePointsKey(const int i, const int j) const override
const LibUtilities::BasisKey v_GetTraceBasisKey(const int i, const int k, bool UseGLL=false) const override
NekDouble v_PhysEvalFirstDeriv(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
int v_GetVertexMap(int localVertexId, bool useCoeffPacking=false) override
int v_CalcNumberOfCoefficients(const std::vector< unsigned int > &nummodes, int &modes_offset) override
int v_GetTraceIntNcoeffs(const int i) const override
StdPrismExp(const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
int v_NumBndryCoeffs() const override
void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray) override
NekDouble v_PhysEvaluateBasis(const Array< OneD, const NekDouble > &coords, int mode) final
void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
LibUtilities::ShapeType v_DetShapeType() const override
Return Shape of region, using ShapeType enum list; i.e. prism.
void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Calculate the inner product of inarray with respect to the basis B=base0*base1*base2 and put into out...
int v_GetNverts() const override
int v_NumDGBndryCoeffs() const override
void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2) override
Calculate the derivative of the physical points.
void v_GetTraceNumModes(const int fid, int &numModes0, int &numModes1, Orientation faceOrient=eDir1FwdDir1_Dir2FwdDir2) override
void v_GetCoords(Array< OneD, NekDouble > &xi_x, Array< OneD, NekDouble > &xi_y, Array< OneD, NekDouble > &xi_z) override
void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray) override
DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey) override
void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
int v_GetNtraces() const 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
int GetMode(int I, int J, int K)
Compute the local mode number in the expansion for a particular tensorial combination.
int v_GetTraceNcoeffs(const int i) 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
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_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_GetEdgeNcoeffs(const int i) const override
void v_LocCollapsedToLocCoord(const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi) override
void v_GetInteriorMap(Array< OneD, unsigned int > &outarray) override
void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) 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_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override
int v_GetTraceNumPoints(const int i) const override
bool v_IsBoundaryInteriorExpansion() const override
int v_GetNedges() const override
void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Inner product of inarray over region with respect to the object's default expansion basis; output in ...
void v_LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta) override
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
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)
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.
@ eModified_B
Principle Modified Functions .
@ P
Monomial polynomials .
@ eOrtho_A
Principle Orthogonal Functions .
@ eModified_C
Principle Modified Functions .
@ eGLL_Lagrange
Lagrange for SEM basis .
@ eOrtho_C
Principle Orthogonal Functions .
@ eOrtho_B
Principle Orthogonal Functions .
@ eModified_A
Principle Modified Functions .
static const NekDouble kNekZeroTol
std::shared_ptr< StdPrismExp > StdPrismExpSharedPtr
LibUtilities::BasisKey EvaluateQuadFaceBasisKey(const int facedir, const LibUtilities::BasisSharedPtr &faceDirBasis)
@ eFactorSVVDGKerDiffCoeff
@ eFactorSVVPowerKerDiffCoeff
LibUtilities::BasisKey EvaluateTriFaceBasisKey(const int facedir, const LibUtilities::BasisSharedPtr &faceDirBasis, bool UseGLL)
const int kSVVDGFiltermodesmin
const int kSVVDGFiltermodesmax
const NekDouble kSVVDGFilter[9][11]
@ ePhysInterpToEquiSpaced
@ eDir1BwdDir2_Dir2BwdDir1
@ eDir1BwdDir1_Dir2BwdDir2
@ eDir1BwdDir2_Dir2FwdDir1
@ eDir1FwdDir1_Dir2BwdDir2
@ eDir1BwdDir1_Dir2FwdDir2
@ eDir1FwdDir2_Dir2FwdDir1
@ eDir1FwdDir2_Dir2BwdDir1
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 Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Svtvp (scalar times vector plus vector): z = alpha*x + y.
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)
scalarT< T > max(scalarT< T > lhs, scalarT< T > rhs)
scalarT< T > min(scalarT< T > lhs, scalarT< T > rhs)
scalarT< T > sqrt(scalarT< T > in)