35 #include <boost/core/ignore_unused.hpp>
47 StdPrismExp::StdPrismExp()
55 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
58 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
62 "order in 'a' direction is higher than order in 'c' direction");
98 int Qx =
m_base[0]->GetNumPoints();
99 int Qy =
m_base[1]->GetNumPoints();
100 int Qz =
m_base[2]->GetNumPoints();
101 int Qtot = Qx * Qy * Qz;
106 eta_x =
m_base[0]->GetZ();
107 eta_z =
m_base[2]->GetZ();
111 bool Do_1 = (out_dxi1.size() > 0) ?
true :
false;
112 bool Do_3 = (out_dxi3.size() > 0) ?
true :
false;
131 for (k = 0; k < Qz; ++k)
134 &dEta_bar1[0] + k * Qx * Qy, 1,
135 &out_dxi1[0] + k * Qx * Qy, 1);
142 for (k = 0; k < Qz; ++k)
145 &dEta_bar1[0] + k * Qx * Qy, 1,
146 &dEta_bar1[0] + k * Qx * Qy, 1);
150 for (i = 0; i < Qx; ++i)
152 Vmath::Svtvp(Qz * Qy, 1.0 + eta_x[i], &dEta_bar1[0] + i, Qx,
153 &out_dxi3[0] + i, Qx, &out_dxi3[0] + i, Qx);
187 ASSERTL1(
false,
"input dir is out of range");
243 "Basis[1] is not a general tensor type");
247 "Basis[2] is not a general tensor type");
249 if (
m_base[0]->Collocation() &&
m_base[1]->Collocation() &&
254 inarray, 1, outarray, 1);
265 int nquad1 =
m_base[1]->GetNumPoints();
266 int nquad2 =
m_base[2]->GetNumPoints();
267 int order0 =
m_base[0]->GetNumModes();
268 int order1 =
m_base[1]->GetNumModes();
271 nquad1 * nquad2 * order0);
274 m_base[2]->GetBdata(), inarray, outarray, wsp,
true,
284 bool doCheckCollDir0,
bool doCheckCollDir1,
bool doCheckCollDir2)
286 boost::ignore_unused(doCheckCollDir0, doCheckCollDir1, doCheckCollDir2);
289 int nquad0 =
m_base[0]->GetNumPoints();
290 int nquad1 =
m_base[1]->GetNumPoints();
291 int nquad2 =
m_base[2]->GetNumPoints();
292 int nummodes0 =
m_base[0]->GetNumModes();
293 int nummodes1 =
m_base[1]->GetNumModes();
294 int nummodes2 =
m_base[2]->GetNumModes();
298 for (i = mode = 0; i < nummodes0; ++i)
300 Blas::Dgemm(
'N',
'N', nquad2, nummodes1, nummodes2 - i, 1.0,
301 base2.get() + mode * nquad2, nquad2,
302 inarray.get() + mode * nummodes1, nummodes2 - i, 0.0,
303 tmp0.get() + i * nquad2 * nummodes1, nquad2);
304 mode += nummodes2 - i;
309 for (i = 0; i < nummodes1; i++)
312 base2.get() + nquad2, 1,
313 tmp0.get() + nquad2 * (nummodes1 + i), 1);
317 for (i = 0; i < nummodes0; i++)
319 Blas::Dgemm(
'N',
'T', nquad1, nquad2, nummodes1, 1.0, base1.get(),
320 nquad1, tmp0.get() + i * nquad2 * nummodes1, nquad2, 0.0,
321 tmp1.get() + i * nquad2 * nquad1, nquad1);
324 Blas::Dgemm(
'N',
'T', nquad0, nquad2 * nquad1, nummodes0, 1.0, base0.get(),
325 nquad0, tmp1.get(), nquad2 * nquad1, 0.0, outarray.get(),
353 out = (*matsys) * in;
391 "Basis[1] is not a general tensor type");
395 "Basis[2] is not a general tensor type");
397 if (
m_base[0]->Collocation() &&
m_base[1]->Collocation())
411 int nquad1 =
m_base[1]->GetNumPoints();
412 int nquad2 =
m_base[2]->GetNumPoints();
413 int order0 =
m_base[0]->GetNumModes();
414 int order1 =
m_base[1]->GetNumModes();
418 if (multiplybyweights)
425 tmp, outarray, wsp,
true,
true,
true);
431 inarray, outarray, wsp,
true,
true,
true);
441 bool doCheckCollDir0,
bool doCheckCollDir1,
bool doCheckCollDir2)
443 boost::ignore_unused(doCheckCollDir0, doCheckCollDir1, doCheckCollDir2);
447 const int nquad0 =
m_base[0]->GetNumPoints();
448 const int nquad1 =
m_base[1]->GetNumPoints();
449 const int nquad2 =
m_base[2]->GetNumPoints();
450 const int order0 =
m_base[0]->GetNumModes();
451 const int order1 =
m_base[1]->GetNumModes();
452 const int order2 =
m_base[2]->GetNumModes();
456 ASSERTL1(wsp.size() >= nquad1 * nquad2 * order0 + nquad2 * order0 * order1,
457 "Insufficient workspace size");
463 Blas::Dgemm(
'T',
'N', nquad1 * nquad2, order0, nquad0, 1.0, inarray.get(),
464 nquad0, base0.get(), nquad0, 0.0, tmp0.get(), nquad1 * nquad2);
467 Blas::Dgemm(
'T',
'N', nquad2 * order0, order1, nquad1, 1.0, tmp0.get(),
468 nquad1, base1.get(), nquad1, 0.0, tmp1.get(), nquad2 * order0);
471 for (mode = i = 0; i < order0; ++i)
473 Blas::Dgemm(
'T',
'N', order2 - i, order1, nquad2, 1.0,
474 base2.get() + mode * nquad2, nquad2,
475 tmp1.get() + i * nquad2, nquad2 * order0, 0.0,
476 outarray.get() + mode * order1, order2 - i);
484 for (i = 0; i < order1; ++i)
489 tmp1.get() + i * order0 * nquad2 + nquad2, 1);
509 ASSERTL0(dir >= 0 && dir <= 2,
"input dir is out of range");
512 int order0 =
m_base[0]->GetNumModes();
513 int order1 =
m_base[1]->GetNumModes();
514 int nquad0 =
m_base[0]->GetNumPoints();
515 int nquad1 =
m_base[1]->GetNumPoints();
516 int nquad2 =
m_base[2]->GetNumPoints();
526 for (i = 0; i < nquad0; ++i)
528 gfac0[i] = 0.5 * (1 + z0[i]);
532 for (i = 0; i < nquad2; ++i)
534 gfac2[i] = 2.0 / (1 - z2[i]);
540 for (i = 0; i < nquad2; ++i)
543 &inarray[0] + i * nquad0 * nquad1, 1,
544 &tmp0[0] + i * nquad0 * nquad1, 1);
555 m_base[2]->GetBdata(), tmp0, outarray, wsp,
true,
true,
true);
563 m_base[2]->GetBdata(), tmp0, outarray, wsp,
true,
true,
true);
572 for (i = 0; i < nquad1 * nquad2; ++i)
574 Vmath::Vmul(nquad0, &gfac0[0], 1, &tmp0[0] + i * nquad0, 1,
575 &tmp0[0] + i * nquad0, 1);
580 m_base[2]->GetBdata(), tmp0, tmp1, wsp,
true,
true,
true);
585 m_base[2]->GetDbdata(), tmp0, outarray, wsp,
true,
true,
true);
615 eta[0] = 2.0 * (1.0 + xi[0]) / d2 - 1.0;
621 xi[0] = (1.0 + eta[0]) * (1.0 - eta[2]) * 0.5 - 1.0;
638 for (
int k = 0; k < Qz; ++k)
640 for (
int j = 0; j < Qy; ++j)
642 for (
int i = 0; i < Qx; ++i)
644 int s = i + Qx * (j + Qy * k);
645 xi_x[s] = (1.0 - eta_z[k]) * (1.0 + etaBar_x[i]) / 2.0 - 1.0;
656 std::array<NekDouble, 3> &firstOrderDerivs)
663 if ((1 - coll[2]) < 1e-5)
667 EphysDeriv2(totPoints);
668 PhysDeriv(inarray, EphysDeriv0, EphysDeriv1, EphysDeriv2);
671 I[0] =
GetBase()[0]->GetI(coll);
672 I[1] =
GetBase()[1]->GetI(coll + 1);
673 I[2] =
GetBase()[2]->GetI(coll + 2);
683 NekDouble dEta_bar1 = firstOrderDerivs[0];
686 firstOrderDerivs[0] = fac * dEta_bar1;
689 fac = 1.0 / (1.0 - coll[2]);
690 dEta_bar1 = fac * dEta_bar1;
694 firstOrderDerivs[2] += fac * dEta_bar1;
712 const int nm1 =
m_base[1]->GetNumModes();
713 const int nm2 =
m_base[2]->GetNumModes();
714 const int b = 2 * nm2 + 1;
716 const int mode0 = floor(0.5 * (b -
sqrt(b * b - 8.0 * mode / nm1)));
718 mode - nm1 * (mode0 * (nm2 - 1) + 1 - (mode0 - 2) * (mode0 - 1) / 2);
719 const int mode1 = tmp / (nm2 - mode0);
720 const int mode2 = tmp % (nm2 - mode0);
722 if (mode0 == 0 && mode2 == 1 &&
726 return StdExpansion::BaryEvaluateBasis<1>(coll[1], mode1) *
727 StdExpansion::BaryEvaluateBasis<2>(coll[2], mode2);
731 return StdExpansion::BaryEvaluateBasis<0>(coll[0], mode0) *
732 StdExpansion::BaryEvaluateBasis<1>(coll[1], mode1) *
733 StdExpansion::BaryEvaluateBasis<2>(coll[2], mode2);
740 int nummodes[3] = {
m_base[0]->GetNumModes(),
m_base[1]->GetNumModes(),
741 m_base[2]->GetNumModes()};
747 numModes0 = nummodes[0];
748 numModes1 = nummodes[1];
755 numModes0 = nummodes[1];
756 numModes1 = nummodes[2];
763 numModes0 = nummodes[0];
764 numModes1 = nummodes[2];
771 std::swap(numModes0, numModes1);
777 ASSERTL2(i >= 0 && i <= 8,
"edge id is out of range");
779 if (i == 0 || i == 2)
783 else if (i == 1 || i == 3 || i == 8)
825 "BasisType is not a boundary interior form");
828 "BasisType is not a boundary interior form");
831 "BasisType is not a boundary interior form");
833 int P =
m_base[0]->GetNumModes();
834 int Q =
m_base[1]->GetNumModes();
835 int R =
m_base[2]->GetNumModes();
844 "BasisType is not a boundary interior form");
847 "BasisType is not a boundary interior form");
850 "BasisType is not a boundary interior form");
852 int P =
m_base[0]->GetNumModes() - 1;
853 int Q =
m_base[1]->GetNumModes() - 1;
854 int R =
m_base[2]->GetNumModes() - 1;
856 return (
P + 1) * (Q + 1)
857 + 2 * (Q + 1) * (R + 1)
858 + 2 * (R + 1) +
P * (1 + 2 * R -
P);
863 ASSERTL2(i >= 0 && i <= 4,
"face id is out of range");
868 else if (i == 1 || i == 3)
871 return Q + 1 + (
P * (1 + 2 * Q -
P)) / 2;
881 ASSERTL2(i >= 0 && i <= 4,
"face id is out of range");
891 else if (i == 1 || i == 3)
893 return Pi * (2 * Ri - Pi - 1) / 2;
903 ASSERTL2(i >= 0 && i <= 4,
"face id is out of range");
907 return m_base[0]->GetNumPoints() *
m_base[1]->GetNumPoints();
909 else if (i == 1 || i == 3)
911 return m_base[0]->GetNumPoints() *
m_base[2]->GetNumPoints();
915 return m_base[1]->GetNumPoints() *
m_base[2]->GetNumPoints();
922 ASSERTL2(i >= 0 && i <= 4,
"face id is out of range");
923 ASSERTL2(j == 0 || j == 1,
"face direction is out of range");
927 return m_base[j]->GetPointsKey();
929 else if (i == 1 || i == 3)
931 return m_base[2 * j]->GetPointsKey();
935 return m_base[j + 1]->GetPointsKey();
942 ASSERTL2(i >= 0 && i <= 4,
"face id is out of range");
943 ASSERTL2(k >= 0 && k <= 1,
"basis key id is out of range");
951 m_base[k]->GetNumModes());
958 m_base[k + 1]->GetNumModes());
965 m_base[2 * k]->GetNumModes());
975 const std::vector<unsigned int> &nummodes,
int &modes_offset)
978 nummodes[modes_offset], nummodes[modes_offset + 1],
979 nummodes[modes_offset + 2]);
1001 "Mapping not defined for this type of basis");
1005 if (useCoeffPacking ==
true)
1028 ASSERTL0(
false,
"local vertex id must be between 0 and 5");
1054 ASSERTL0(
false,
"local vertex id must be between 0 and 5");
1065 "BasisType is not a boundary interior form");
1068 "BasisType is not a boundary interior form");
1071 "BasisType is not a boundary interior form");
1073 int P =
m_base[0]->GetNumModes() - 1,
p;
1074 int Q =
m_base[1]->GetNumModes() - 1, q;
1075 int R =
m_base[2]->GetNumModes() - 1, r;
1079 if (outarray.size() != nIntCoeffs)
1087 for (
p = 2;
p <=
P; ++
p)
1089 for (q = 2; q <= Q; ++q)
1091 for (r = 1; r <= R -
p; ++r)
1093 outarray[idx++] =
GetMode(
p, q, r);
1103 "BasisType is not a boundary interior form");
1106 "BasisType is not a boundary interior form");
1109 "BasisType is not a boundary interior form");
1111 int P =
m_base[0]->GetNumModes() - 1,
p;
1112 int Q =
m_base[1]->GetNumModes() - 1, q;
1113 int R =
m_base[2]->GetNumModes() - 1, r;
1118 if (maparray.size() != nBnd)
1124 for (
p = 0;
p <=
P; ++
p)
1129 for (q = 0; q <= Q; ++q)
1131 for (r = 0; r <= R -
p; ++r)
1133 maparray[idx++] =
GetMode(
p, q, r);
1141 for (q = 0; q <= Q; ++q)
1145 for (r = 0; r <= R -
p; ++r)
1147 maparray[idx++] =
GetMode(
p, q, r);
1152 maparray[idx++] =
GetMode(
p, q, 0);
1163 "Method only implemented if BasisType is identical"
1164 "in x and y directions");
1167 "Method only implemented for Modified_A BasisType"
1168 "(x and y direction) and Modified_B BasisType (z "
1170 int p, q, r, idx = 0;
1177 Q =
m_base[1]->GetNumModes();
1182 Q =
m_base[2]->GetNumModes();
1187 Q =
m_base[2]->GetNumModes();
1190 ASSERTL0(
false,
"fid must be between 0 and 4");
1193 if (maparray.size() !=
P * Q)
1203 for (q = 0; q < Q; ++q)
1205 for (
p = 0;
p <
P; ++
p)
1212 for (
p = 0;
p <
P; ++
p)
1214 for (r = 0; r < Q -
p; ++r)
1216 maparray[idx++] =
GetMode(
p, 0, r);
1221 for (q = 0; q <
P; ++q)
1223 maparray[q] =
GetMode(1, q, 0);
1225 for (q = 0; q <
P; ++q)
1227 maparray[
P + q] =
GetMode(0, q, 1);
1229 for (r = 1; r < Q - 1; ++r)
1231 for (q = 0; q <
P; ++q)
1233 maparray[(r + 1) *
P + q] =
GetMode(1, q, r);
1238 for (
p = 0;
p <
P; ++
p)
1240 for (r = 0; r < Q -
p; ++r)
1242 maparray[idx++] =
GetMode(
p, 1, r);
1247 for (r = 0; r < Q; ++r)
1249 for (q = 0; q <
P; ++q)
1251 maparray[r *
P + q] =
GetMode(0, q, r);
1256 ASSERTL0(
false,
"Face to element map unavailable.");
1266 "Method only implemented if BasisType is identical"
1267 "in x and y directions");
1270 "Method only implemented for Modified_A BasisType"
1271 "(x and y direction) and Modified_B BasisType (z "
1274 int i, j, k,
p, r, nFaceCoeffs, idx = 0;
1275 int nummodesA = 0, nummodesB = 0;
1280 nummodesA =
m_base[0]->GetNumModes();
1281 nummodesB =
m_base[1]->GetNumModes();
1285 nummodesA =
m_base[0]->GetNumModes();
1286 nummodesB =
m_base[2]->GetNumModes();
1290 nummodesA =
m_base[1]->GetNumModes();
1291 nummodesB =
m_base[2]->GetNumModes();
1294 ASSERTL0(
false,
"fid must be between 0 and 4");
1303 else if (fid == 1 || fid == 3)
1305 nFaceCoeffs =
P * (2 * Q -
P + 1) / 2;
1309 nFaceCoeffs =
P * Q;
1313 if (maparray.size() != nFaceCoeffs)
1318 if (signarray.size() != nFaceCoeffs)
1324 fill(signarray.get(), signarray.get() + nFaceCoeffs, 1);
1327 int minPA = min(nummodesA,
P);
1328 int minQB = min(nummodesB, Q);
1330 if (fid == 1 || fid == 3)
1337 for (j = 0; j < minPA; ++j)
1340 for (k = 0; k < minQB - j; ++k, ++cnt)
1342 maparray[idx++] = cnt;
1345 cnt += nummodesB - minQB;
1348 for (k = nummodesB - j; k < Q - j; ++k)
1350 signarray[idx] = 0.0;
1351 maparray[idx++] = maparray[0];
1355 for (j = minPA; j < nummodesA; ++j)
1358 for (k = 0; k < minQB-j; ++k, ++cnt)
1360 maparray[idx++] = cnt;
1363 cnt += nummodesB-minQB;
1366 for (k = nummodesB-j; k < Q-j; ++k)
1368 signarray[idx] = 0.0;
1369 maparray[idx++] = maparray[0];
1373 for (j = nummodesA; j <
P; ++j)
1375 for (k = 0; k < Q - j; ++k)
1377 signarray[idx] = 0.0;
1378 maparray[idx++] = maparray[0];
1387 for (
p = 0;
p <
P; ++
p)
1389 for (r = 0; r < Q -
p; ++r, idx++)
1393 signarray[idx] =
p % 2 ? -1 : 1;
1398 swap(maparray[0], maparray[Q]);
1399 for (i = 1; i < Q - 1; ++i)
1401 swap(maparray[i + 1], maparray[Q + i]);
1411 for (i = 0; i < Q; i++)
1413 for (j = 0; j <
P; j++)
1417 arrayindx[i *
P + j] = i *
P + j;
1421 arrayindx[i *
P + j] = j * Q + i;
1428 for (j = 0; j <
P; ++j)
1431 for (k = 0; k < Q; k++)
1433 maparray[arrayindx[j + k *
P]] = j + k * nummodesA;
1436 for (k = nummodesB; k < Q; ++k)
1438 signarray[arrayindx[j + k *
P]] = 0.0;
1439 maparray[arrayindx[j + k *
P]] = maparray[0];
1443 for (j = nummodesA; j <
P; ++j)
1445 for (k = 0; k < Q; ++k)
1447 signarray[arrayindx[j + k *
P]] = 0.0;
1448 maparray[arrayindx[j + k *
P]] = maparray[0];
1463 for (i = 3; i < Q; i += 2)
1465 for (j = 0; j <
P; j++)
1467 signarray[arrayindx[i *
P + j]] *= -1;
1471 for (i = 0; i <
P; i++)
1473 swap(maparray[i], maparray[i +
P]);
1474 swap(signarray[i], signarray[i +
P]);
1479 for (i = 0; i < Q; i++)
1481 for (j = 3; j <
P; j += 2)
1483 signarray[arrayindx[i *
P + j]] *= -1;
1487 for (i = 0; i < Q; i++)
1489 swap(maparray[i], maparray[i + Q]);
1490 swap(signarray[i], signarray[i + Q]);
1502 for (i = 0; i < Q; i++)
1504 for (j = 3; j <
P; j += 2)
1506 signarray[arrayindx[i *
P + j]] *= -1;
1510 for (i = 0; i < Q; i++)
1512 swap(maparray[i *
P], maparray[i *
P + 1]);
1513 swap(signarray[i *
P], signarray[i *
P + 1]);
1518 for (i = 3; i < Q; i += 2)
1520 for (j = 0; j <
P; j++)
1522 signarray[arrayindx[i *
P + j]] *= -1;
1526 for (i = 0; i <
P; i++)
1528 swap(maparray[i * Q], maparray[i * Q + 1]);
1529 swap(signarray[i * Q], signarray[i * Q + 1]);
1542 const int P =
m_base[0]->GetNumModes() - 1;
1543 const int Q =
m_base[1]->GetNumModes() - 1;
1544 const int R =
m_base[2]->GetNumModes() - 1;
1547 if (maparray.size() != nEdgeIntCoeffs)
1552 if (signarray.size() != nEdgeIntCoeffs)
1558 fill(signarray.get(), signarray.get() + nEdgeIntCoeffs, 1);
1568 for (i = 2; i <=
P; ++i)
1570 maparray[i - 2] =
GetMode(i, 0, 0);
1575 for (i = 2; i <= Q; ++i)
1577 maparray[i - 2] =
GetMode(1, i, 0);
1584 for (i = 2; i <=
P; ++i)
1586 maparray[i - 2] =
GetMode(i, 1, 0);
1593 for (i = 2; i <= Q; ++i)
1595 maparray[i - 2] =
GetMode(0, i, 0);
1600 for (i = 2; i <= R; ++i)
1602 maparray[i - 2] =
GetMode(0, 0, i);
1607 for (i = 1; i <= R - 1; ++i)
1609 maparray[i - 1] =
GetMode(1, 0, i);
1614 for (i = 1; i <= R - 1; ++i)
1616 maparray[i - 1] =
GetMode(1, 1, i);
1621 for (i = 2; i <= R; ++i)
1623 maparray[i - 2] =
GetMode(0, 1, i);
1628 for (i = 2; i <= Q; ++i)
1630 maparray[i - 2] =
GetMode(0, i, 1);
1635 ASSERTL0(
false,
"Edge not defined.");
1641 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1652 const int P =
m_base[0]->GetNumModes() - 1;
1653 const int Q =
m_base[1]->GetNumModes() - 1;
1654 const int R =
m_base[2]->GetNumModes() - 1;
1656 int p, q, r, idx = 0;
1662 if (maparray.size() != nFaceIntCoeffs)
1667 if (signarray.size() != nFaceIntCoeffs)
1673 fill(signarray.get(), signarray.get() + nFaceIntCoeffs, 1);
1679 if (fid != 1 && fid != 3)
1692 for (i = 0; i < nummodesB; i++)
1694 for (j = 0; j < nummodesA; j++)
1698 arrayindx[i * nummodesA + j] = i * nummodesA + j;
1702 arrayindx[i * nummodesA + j] = j * nummodesB + i;
1711 for (q = 2; q <= Q; ++q)
1713 for (
p = 2;
p <=
P; ++
p)
1715 maparray[arrayindx[(q - 2) * nummodesA + (
p - 2)]] =
1722 for (
p = 2;
p <=
P; ++
p)
1724 for (r = 1; r <= R -
p; ++r)
1728 signarray[idx] =
p % 2 ? -1 : 1;
1730 maparray[idx++] =
GetMode(
p, 0, r);
1736 for (r = 1; r <= R - 1; ++r)
1738 for (q = 2; q <= Q; ++q)
1740 maparray[arrayindx[(r - 1) * nummodesA + (q - 2)]] =
1747 for (
p = 2;
p <=
P; ++
p)
1749 for (r = 1; r <= R -
p; ++r)
1753 signarray[idx] =
p % 2 ? -1 : 1;
1755 maparray[idx++] =
GetMode(
p, 1, r);
1761 for (r = 2; r <= R; ++r)
1763 for (q = 2; q <= Q; ++q)
1765 maparray[arrayindx[(r - 2) * nummodesA + (q - 2)]] =
1772 ASSERTL0(
false,
"Face interior map not available.");
1777 if (fid == 1 || fid == 3)
1787 for (i = 1; i < nummodesB; i += 2)
1789 for (j = 0; j < nummodesA; j++)
1791 signarray[arrayindx[i * nummodesA + j]] *= -1;
1797 for (i = 0; i < nummodesB; i++)
1799 for (j = 1; j < nummodesA; j += 2)
1801 signarray[arrayindx[i * nummodesA + j]] *= -1;
1814 for (i = 0; i < nummodesB; i++)
1816 for (j = 1; j < nummodesA; j += 2)
1818 signarray[arrayindx[i * nummodesA + j]] *= -1;
1824 for (i = 1; i < nummodesB; i += 2)
1826 for (j = 0; j < nummodesA; j++)
1828 signarray[arrayindx[i * nummodesA + j]] *= -1;
1850 int nq0 =
m_base[0]->GetNumPoints();
1851 int nq1 =
m_base[1]->GetNumPoints();
1852 int nq2 =
m_base[2]->GetNumPoints();
1862 nq = max(nq0, max(nq1, nq2));
1875 for (
int i = 0; i < nq; ++i)
1877 for (
int j = 0; j < nq; ++j)
1879 for (
int k = 0; k < nq - i; ++k, ++cnt)
1882 coords[cnt][0] = -1.0 + 2 * k / (
NekDouble)(nq - 1);
1883 coords[cnt][1] = -1.0 + 2 * j / (
NekDouble)(nq - 1);
1884 coords[cnt][2] = -1.0 + 2 * i / (
NekDouble)(nq - 1);
1889 for (
int i = 0; i < neq; ++i)
1893 I[0] =
m_base[0]->GetI(coll);
1894 I[1] =
m_base[1]->GetI(coll + 1);
1895 I[2] =
m_base[2]->GetI(coll + 2);
1899 for (
int k = 0; k < nq2; ++k)
1901 for (
int j = 0; j < nq1; ++j)
1904 fac = (I[1]->GetPtr())[j] * (I[2]->GetPtr())[k];
1908 Mat->GetRawPtr() + k * nq0 * nq1 * neq +
1951 int Q =
m_base[1]->GetNumModes() - 1;
1952 int R =
m_base[2]->GetNumModes() - 1;
1956 (Q + 1) * (
p * R + 1 -
1957 (
p - 2) * (
p - 1) / 2);
1965 int nquad0 =
m_base[0]->GetNumPoints();
1966 int nquad1 =
m_base[1]->GetNumPoints();
1967 int nquad2 =
m_base[2]->GetNumPoints();
1976 for (i = 0; i < nquad1 * nquad2; ++i)
1978 Vmath::Vmul(nquad0, inarray.get() + i * nquad0, 1, w0.get(), 1,
1979 outarray.get() + i * nquad0, 1);
1983 for (j = 0; j < nquad2; ++j)
1985 for (i = 0; i < nquad1; ++i)
1988 &outarray[0] + i * nquad0 + j * nquad0 * nquad1, 1);
1998 case LibUtilities::eGaussRadauMAlpha1Beta0:
1999 for (i = 0; i < nquad2; ++i)
2002 &outarray[0] + i * nquad0 * nquad1, 1);
2007 for (i = 0; i < nquad2; ++i)
2009 Blas::Dscal(nquad0 * nquad1, 0.5 * (1 - z2[i]) * w2[i],
2010 &outarray[0] + i * nquad0 * nquad1, 1);
2020 int qa =
m_base[0]->GetNumPoints();
2021 int qb =
m_base[1]->GetNumPoints();
2022 int qc =
m_base[2]->GetNumPoints();
2023 int nmodes_a =
m_base[0]->GetNumModes();
2024 int nmodes_b =
m_base[1]->GetNumModes();
2025 int nmodes_c =
m_base[2]->GetNumModes();
2037 int i, j, k, cnt = 0;
2040 OrthoExp.
FwdTrans(array, orthocoeffs);
2050 for (
int i = 0; i < nmodes_a; ++i)
2052 for (
int j = 0; j < nmodes_b; ++j)
2055 pow((1.0 * i) / (nmodes_a - 1), cutoff * nmodes_a),
2056 pow((1.0 * j) / (nmodes_b - 1), cutoff * nmodes_b));
2058 for (
int k = 0; k < nmodes_c - i; ++k)
2061 std::max(fac1, pow((1.0 * k) / (nmodes_c - 1),
2062 cutoff * nmodes_c));
2064 orthocoeffs[cnt] *= SvvDiffCoeff * fac;
2080 max_abc = max(max_abc, 0);
2083 for (
int i = 0; i < nmodes_a; ++i)
2085 for (
int j = 0; j < nmodes_b; ++j)
2087 int maxij = max(i, j);
2089 for (
int k = 0; k < nmodes_c - i; ++k)
2091 int maxijk = max(maxij, k);
2113 int cutoff_a = (int)(SVVCutOff * nmodes_a);
2114 int cutoff_b = (int)(SVVCutOff * nmodes_b);
2115 int cutoff_c = (int)(SVVCutOff * nmodes_c);
2119 int nmodes = min(min(nmodes_a, nmodes_b), nmodes_c);
2120 NekDouble cutoff = min(min(cutoff_a, cutoff_b), cutoff_c);
2123 for (i = 0; i < nmodes_a; ++i)
2125 for (j = 0; j < nmodes_b; ++j)
2127 for (k = 0; k < nmodes_c - i; ++k)
2129 if (j >= cutoff || i + k >= cutoff)
2133 exp(-(i + k - nmodes) * (i + k - nmodes) /
2134 ((
NekDouble)((i + k - cutoff + epsilon) *
2135 (i + k - cutoff + epsilon)))) *
2136 exp(-(j - nmodes) * (j - nmodes) /
2138 (j - cutoff + epsilon)))));
2142 orthocoeffs[cnt] *= 0.0;
2151 OrthoExp.
BwdTrans(orthocoeffs, array);
2158 int nquad0 =
m_base[0]->GetNumPoints();
2159 int nquad1 =
m_base[1]->GetNumPoints();
2160 int nquad2 =
m_base[2]->GetNumPoints();
2161 int nqtot = nquad0 * nquad1 * nquad2;
2162 int nmodes0 =
m_base[0]->GetNumModes();
2163 int nmodes1 =
m_base[1]->GetNumModes();
2164 int nmodes2 =
m_base[2]->GetNumModes();
2165 int numMax = nmodes0;
2186 bortho0, bortho1, bortho2);
2189 OrthoPrismExp->FwdTrans(phys_tmp, coeff);
2192 for (u = 0; u < numMin; ++u)
2194 for (i = 0; i < numMin; ++i)
2197 tmp2 = coeff_tmp1 + cnt, 1);
2201 for (i = numMin; i < numMax; ++i)
2207 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
virtual 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
virtual void v_GetEdgeInteriorToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2) override
virtual LibUtilities::PointsKey v_GetTracePointsKey(const int i, const int j) const override
virtual int v_GetVertexMap(int localVertexId, bool useCoeffPacking=false) override
virtual int v_CalcNumberOfCoefficients(const std::vector< unsigned int > &nummodes, int &modes_offset) override
virtual int v_GetTraceIntNcoeffs(const int i) const override
virtual int v_NumBndryCoeffs() const override
virtual void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray) override
void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
virtual const LibUtilities::BasisKey v_GetTraceBasisKey(const int i, const int k) const override
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
virtual 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...
virtual int v_GetNverts() const override
virtual 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.
virtual void v_GetTraceNumModes(const int fid, int &numModes0, int &numModes1, Orientation faceOrient=eDir1FwdDir1_Dir2FwdDir2) override
virtual void v_GetCoords(Array< OneD, NekDouble > &xi_x, Array< OneD, NekDouble > &xi_y, Array< OneD, NekDouble > &xi_z) override
virtual void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray) override
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey) override
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
virtual int v_GetNtraces() const override
virtual ~StdPrismExp() override
virtual 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.
virtual int v_GetTraceNcoeffs(const int i) const override
virtual 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
virtual 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...
NekDouble v_PhysEvaluateBasis(const Array< OneD, const NekDouble > &coords, int mode) final override
virtual NekDouble v_PhysEvaluate(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
virtual 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
virtual int v_GetEdgeNcoeffs(const int i) const override
virtual void v_LocCollapsedToLocCoord(const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi) override
virtual void v_GetInteriorMap(Array< OneD, unsigned int > &outarray) override
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
virtual 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
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override
virtual int v_GetTraceNumPoints(const int i) const override
virtual bool v_IsBoundaryInteriorExpansion() const override
virtual int v_GetNedges() const override
virtual 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 ...
virtual 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)
int getNumberOfCoefficients(int Na)
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.
@ eModified_B
Principle Modified Functions .
@ 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 EvaluateTriFaceBasisKey(const int facedir, const LibUtilities::BasisType faceDirBasisType, const int numpoints, const int nummodes)
@ eFactorSVVDGKerDiffCoeff
@ eFactorSVVPowerKerDiffCoeff
const int kSVVDGFiltermodesmin
const int kSVVDGFiltermodesmax
LibUtilities::BasisKey EvaluateQuadFaceBasisKey(const int facedir, const LibUtilities::BasisType faceDirBasisType, const int numpoints, const int nummodes)
const NekDouble kSVVDGFilter[9][11]
@ ePhysInterpToEquiSpaced
@ eDir1BwdDir2_Dir2BwdDir1
@ eDir1BwdDir1_Dir2BwdDir2
@ eDir1BwdDir2_Dir2FwdDir1
@ eDir1FwdDir1_Dir2BwdDir2
@ eDir1BwdDir1_Dir2FwdDir2
@ eDir1FwdDir2_Dir2FwdDir1
@ eDir1FwdDir2_Dir2BwdDir1
The above copyright notice and this permission notice shall be included.
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 > sqrt(scalarT< T > in)