45 :
StdExpansion(LibUtilities::StdTetData::getNumberOfCoefficients(
46 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
49 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
53 "order in 'a' direction is higher than order "
56 "order in 'a' direction is higher than order "
59 "order in 'b' direction is higher than order "
88 int Q0 =
m_base[0]->GetNumPoints();
89 int Q1 =
m_base[1]->GetNumPoints();
90 int Q2 =
m_base[2]->GetNumPoints();
91 int Qtot = Q0 * Q1 * Q2;
98 bool Do_2 = (out_dxi2.size() > 0) ?
true :
false;
99 bool Do_1 = (out_dxi1.size() > 0) ?
true :
false;
116 eta_0 =
m_base[0]->GetZ();
117 eta_1 =
m_base[1]->GetZ();
118 eta_2 =
m_base[2]->GetZ();
124 for (
int k = 0; k < Q2; ++k)
126 for (
int j = 0; j < Q1; ++j, dEta0 += Q0)
128 Vmath::Smul(Q0, 2.0 / (1.0 - eta_1[j]), dEta0, 1, dEta0, 1);
130 fac = 1.0 / (1.0 - eta_2[k]);
131 Vmath::Smul(Q0 * Q1, fac, &out_dEta0[0] + k * Q0 * Q1, 1,
132 &out_dEta0[0] + k * Q0 * Q1, 1);
135 if (out_dxi0.size() > 0)
147 for (
int k = 0; k < Q1 * Q2; ++k)
149 Vmath::Vmul(Q0, &Fac0[0], 1, &out_dEta0[0] + k * Q0, 1,
150 &out_dEta0[0] + k * Q0, 1);
153 for (
int k = 0; k < Q2; ++k)
156 &out_dEta1[0] + k * Q0 * Q1, 1,
157 &out_dEta1[0] + k * Q0 * Q1, 1);
164 Vmath::Vadd(Qtot, out_dEta0, 1, out_dEta1, 1, out_dxi1, 1);
171 for (
int k = 0; k < Q2; ++k)
173 for (
int j = 0; j < Q1; ++j, dEta1 += Q0)
175 Vmath::Smul(Q0, (1.0 + eta_1[j]) / 2.0, dEta1, 1, dEta1, 1);
182 Vmath::Vadd(Qtot, out_dEta0, 1, out_dEta1, 1, out_dxi2, 1);
183 Vmath::Vadd(Qtot, out_dEta2, 1, out_dxi2, 1, out_dxi2, 1);
220 ASSERTL1(
false,
"input dir is out of range");
264 "Basis[1] is not a general tensor type");
268 "Basis[2] is not a general tensor type");
270 if (
m_base[0]->Collocation() &&
m_base[1]->Collocation() &&
275 inarray, 1, outarray, 1);
289 int nquad1 =
m_base[1]->GetNumPoints();
290 int nquad2 =
m_base[2]->GetNumPoints();
291 int order0 =
m_base[0]->GetNumModes();
292 int order1 =
m_base[1]->GetNumModes();
295 nquad2 * nquad1 * order0);
298 m_base[2]->GetBdata(), inarray, outarray, wsp,
true,
321 [[maybe_unused]]
bool doCheckCollDir0,
322 [[maybe_unused]]
bool doCheckCollDir1,
323 [[maybe_unused]]
bool doCheckCollDir2)
325 int nquad0 =
m_base[0]->GetNumPoints();
326 int nquad1 =
m_base[1]->GetNumPoints();
327 int nquad2 =
m_base[2]->GetNumPoints();
329 int order0 =
m_base[0]->GetNumModes();
330 int order1 =
m_base[1]->GetNumModes();
331 int order2 =
m_base[2]->GetNumModes();
335 tmp + nquad2 * order0 * (2 * order1 - order0 + 1) / 2;
337 int i, j, mode, mode1, cnt;
340 mode = mode1 = cnt = 0;
341 for (i = 0; i < order0; ++i)
343 for (j = 0; j < order1 - i; ++j, ++cnt)
346 base2.data() + mode * nquad2, nquad2,
347 inarray.data() + mode1, 1, 0.0,
348 tmp.data() + cnt * nquad2, 1);
349 mode += order2 - i - j;
350 mode1 += order2 - i - j;
353 for (j = order1 - i; j < order2 - i; ++j)
355 mode += order2 - i - j;
365 Blas::Daxpy(nquad2, inarray[1], base2.data() + nquad2, 1,
366 &tmp[0] + nquad2, 1);
369 Blas::Daxpy(nquad2, inarray[1], base2.data() + nquad2, 1,
370 &tmp[0] + order1 * nquad2, 1);
375 for (i = 0; i < order0; ++i)
377 Blas::Dgemm(
'N',
'T', nquad1, nquad2, order1 - i, 1.0,
378 base1.data() + mode * nquad1, nquad1,
379 tmp.data() + mode * nquad2, nquad2, 0.0,
380 tmp1.data() + i * nquad1 * nquad2, nquad1);
391 for (i = 0; i < nquad2; ++i)
393 Blas::Daxpy(nquad1, tmp[nquad2 + i], base1.data() + nquad1, 1,
394 &tmp1[nquad1 * nquad2] + i * nquad1, 1);
399 Blas::Dgemm(
'N',
'T', nquad0, nquad1 * nquad2, order0, 1.0, base0.data(),
400 nquad0, tmp1.data(), nquad1 * nquad2, 0.0, outarray.data(),
422 out = (*matsys) * in;
464 "Basis[1] is not a general tensor type");
468 "Basis[2] is not a general tensor type");
470 if (
m_base[0]->Collocation() &&
m_base[1]->Collocation())
490 int nquad0 =
m_base[0]->GetNumPoints();
491 int nquad1 =
m_base[1]->GetNumPoints();
492 int nquad2 =
m_base[2]->GetNumPoints();
493 int order0 =
m_base[0]->GetNumModes();
494 int order1 =
m_base[1]->GetNumModes();
497 nquad2 * order0 * (2 * order1 - order0 + 1) / 2);
499 if (multiplybyweights)
506 tmp, outarray, wsp,
true,
true,
true);
512 inarray, outarray, wsp,
true,
true,
true);
522 [[maybe_unused]]
bool doCheckCollDir0,
523 [[maybe_unused]]
bool doCheckCollDir1,
524 [[maybe_unused]]
bool doCheckCollDir2)
526 int nquad0 =
m_base[0]->GetNumPoints();
527 int nquad1 =
m_base[1]->GetNumPoints();
528 int nquad2 =
m_base[2]->GetNumPoints();
530 int order0 =
m_base[0]->GetNumModes();
531 int order1 =
m_base[1]->GetNumModes();
532 int order2 =
m_base[2]->GetNumModes();
537 int i, j, mode, mode1, cnt;
540 Blas::Dgemm(
'T',
'N', nquad1 * nquad2, order0, nquad0, 1.0, inarray.data(),
541 nquad0, base0.data(), nquad0, 0.0, tmp1.data(),
545 for (mode = i = 0; i < order0; ++i)
547 Blas::Dgemm(
'T',
'N', nquad2, order1 - i, nquad1, 1.0,
548 tmp1.data() + i * nquad1 * nquad2, nquad1,
549 base1.data() + mode * nquad1, nquad1, 0.0,
550 tmp2.data() + mode * nquad2, nquad2);
559 Blas::Dgemv(
'T', nquad1, nquad2, 1.0, tmp1.data() + nquad1 * nquad2,
560 nquad1, base1.data() + nquad1, 1, 1.0, tmp2.data() + nquad2,
565 mode = mode1 = cnt = 0;
566 for (i = 0; i < order0; ++i)
568 for (j = 0; j < order1 - i; ++j, ++cnt)
571 base2.data() + mode * nquad2, nquad2,
572 tmp2.data() + cnt * nquad2, 1, 0.0,
573 outarray.data() + mode1, 1);
574 mode += order2 - i - j;
575 mode1 += order2 - i - j;
578 for (j = order1 - i; j < order2 - i; ++j)
580 mode += order2 - i - j;
590 Blas::Ddot(nquad2, base2.data() + nquad2, 1, &tmp2[nquad2], 1);
593 outarray[1] +=
Blas::Ddot(nquad2, base2.data() + nquad2, 1,
594 &tmp2[nquad2 * order1], 1);
616 int nquad0 =
m_base[0]->GetNumPoints();
617 int nquad1 =
m_base[1]->GetNumPoints();
618 int nquad2 =
m_base[2]->GetNumPoints();
619 int nqtot = nquad0 * nquad1 * nquad2;
620 int nmodes0 =
m_base[0]->GetNumModes();
621 int nmodes1 =
m_base[1]->GetNumModes();
622 int wspsize = nquad0 + nquad1 + nquad2 +
max(nqtot,
m_ncoeffs) +
623 nquad1 * nquad2 * nmodes0 +
624 nquad2 * nmodes0 * (2 * nmodes1 - nmodes0 + 1) / 2;
637 for (i = 0; i < nquad0; ++i)
639 gfac0[i] = 0.5 * (1 + z0[i]);
643 for (i = 0; i < nquad1; ++i)
645 gfac1[i] = 2.0 / (1 - z1[i]);
649 for (i = 0; i < nquad2; ++i)
651 gfac2[i] = 2.0 / (1 - z2[i]);
655 for (i = 0; i < nquad1 * nquad2; ++i)
657 Vmath::Smul(nquad0, gfac1[i % nquad1], &inarray[0] + i * nquad0, 1,
658 &tmp0[0] + i * nquad0, 1);
660 for (i = 0; i < nquad2; ++i)
662 Vmath::Smul(nquad0 * nquad1, gfac2[i], &tmp0[0] + i * nquad0 * nquad1,
663 1, &tmp0[0] + i * nquad0 * nquad1, 1);
674 m_base[2]->GetBdata(), tmp0, outarray, wsp,
false,
true,
true);
681 for (i = 0; i < nquad1 * nquad2; ++i)
683 Vmath::Vmul(nquad0, &gfac0[0], 1, &tmp0[0] + i * nquad0, 1,
684 &tmp0[0] + i * nquad0, 1);
689 m_base[2]->GetBdata(), tmp0, tmp3, wsp,
false,
true,
true);
691 for (i = 0; i < nquad2; ++i)
694 &inarray[0] + i * nquad0 * nquad1, 1,
695 &tmp0[0] + i * nquad0 * nquad1, 1);
700 m_base[2]->GetBdata(), tmp0, outarray, wsp,
true,
false,
true);
709 for (i = 0; i < nquad1; ++i)
711 gfac1[i] = (1 + z1[i]) / 2;
714 for (i = 0; i < nquad1 * nquad2; ++i)
716 Vmath::Vmul(nquad0, &gfac0[0], 1, &tmp0[0] + i * nquad0, 1,
717 &tmp0[0] + i * nquad0, 1);
721 m_base[2]->GetBdata(), tmp0, tmp3, wsp,
false,
true,
true);
723 for (i = 0; i < nquad2; ++i)
726 &inarray[0] + i * nquad0 * nquad1, 1,
727 &tmp0[0] + i * nquad0 * nquad1, 1);
729 for (i = 0; i < nquad1 * nquad2; ++i)
731 Vmath::Smul(nquad0, gfac1[i % nquad1], &tmp0[0] + i * nquad0, 1,
732 &tmp0[0] + i * nquad0, 1);
737 m_base[2]->GetBdata(), tmp0, tmp4, wsp,
true,
false,
true);
742 m_base[2]->GetDbdata(), tmp0, outarray, wsp,
true,
true,
false);
752 ASSERTL1(
false,
"input dir is out of range");
789 eta[0] = 2.0 * (1.0 + xi[0]) / d12 - 1.0;
790 eta[1] = 2.0 * (1.0 + xi[1]) / d2 - 1.0;
798 xi[1] = (1.0 + eta[1]) * (1.0 - xi[2]) * 0.5 - 1.0;
799 xi[0] = (1.0 + eta[0]) * (-xi[1] - xi[2]) * 0.5 - 1.0;
815 const int nm1 =
m_base[1]->GetNumModes();
816 const int nm2 =
m_base[2]->GetNumModes();
818 const int b = 2 * nm2 + 1;
819 const int mode0 = floor(0.5 * (b -
sqrt(b * b - 8.0 * mode / nm1)));
821 mode - nm1 * (mode0 * (nm2 - 1) + 1 - (mode0 - 2) * (mode0 - 1) / 2);
822 const int mode1 = tmp / (nm2 - mode0);
823 const int mode2 = tmp % (nm2 - mode0);
832 return StdExpansion::BaryEvaluateBasis<2>(coll[2], 1);
834 else if (mode0 == 0 && mode2 == 1)
836 return StdExpansion::BaryEvaluateBasis<1>(coll[1], 0) *
837 StdExpansion::BaryEvaluateBasis<2>(coll[2], 1);
839 else if (mode0 == 1 && mode1 == 1 && mode2 == 0)
841 return StdExpansion::BaryEvaluateBasis<0>(coll[0], 0) *
842 StdExpansion::BaryEvaluateBasis<1>(coll[1], 1);
846 return StdExpansion::BaryEvaluateBasis<0>(coll[0], mode0) *
847 StdExpansion::BaryEvaluateBasis<1>(coll[1], mode1) *
848 StdExpansion::BaryEvaluateBasis<2>(coll[2], mode2);
854 std::array<NekDouble, 3> &firstOrderDerivs)
861 if ((1 - coll[1]) < 1e-5 || (1 - coll[2]) < 1e-5)
865 EphysDeriv2(totPoints);
866 PhysDeriv(inarray, EphysDeriv0, EphysDeriv1, EphysDeriv2);
869 I[0] =
GetBase()[0]->GetI(coll);
870 I[1] =
GetBase()[1]->GetI(coll + 1);
871 I[2] =
GetBase()[2]->GetI(coll + 2);
879 std::array<NekDouble, 3> interDeriv;
883 NekDouble temp = 2.0 / ((1 - coll[1]) * (1 - coll[2]));
884 interDeriv[0] *= temp;
887 firstOrderDerivs[0] = 2 * interDeriv[0];
894 interDeriv[0] *= fac0;
897 fac0 = 2 / (1 - coll[2]);
898 interDeriv[1] *= fac0;
902 firstOrderDerivs[1] = interDeriv[0] + interDeriv[1];
905 fac0 = (1 + coll[1]) / 2;
906 interDeriv[1] *= fac0;
911 firstOrderDerivs[2] = interDeriv[0] + interDeriv[1] + interDeriv[2];
920 int nummodes[3] = {
m_base[0]->GetNumModes(),
m_base[1]->GetNumModes(),
921 m_base[2]->GetNumModes()};
926 numModes0 = nummodes[0];
927 numModes1 = nummodes[1];
932 numModes0 = nummodes[0];
933 numModes1 = nummodes[2];
939 numModes0 = nummodes[1];
940 numModes1 = nummodes[2];
974 "BasisType is not a boundary interior form");
977 "BasisType is not a boundary interior form");
980 "BasisType is not a boundary interior form");
982 int P =
m_base[0]->GetNumModes();
983 int Q =
m_base[1]->GetNumModes();
984 int R =
m_base[2]->GetNumModes();
993 "BasisType is not a boundary interior form");
996 "BasisType is not a boundary interior form");
999 "BasisType is not a boundary interior form");
1001 int P =
m_base[0]->GetNumModes() - 1;
1002 int Q =
m_base[1]->GetNumModes() - 1;
1003 int R =
m_base[2]->GetNumModes() - 1;
1005 return (Q + 1) +
P * (1 + 2 * Q -
P) / 2
1006 + (R + 1) +
P * (1 + 2 * R -
P) / 2
1007 + 2 * (R + 1) + Q * (1 + 2 * R - Q);
1012 ASSERTL2((i >= 0) && (i <= 3),
"face id is out of range");
1013 int nFaceCoeffs = 0;
1014 int nummodesA, nummodesB,
P, Q;
1020 else if ((i == 1) || (i == 2))
1032 nFaceCoeffs = Q + 1 + (
P * (1 + 2 * Q -
P)) / 2;
1038 ASSERTL2((i >= 0) && (i <= 3),
"face id is out of range");
1039 int Pi =
m_base[0]->GetNumModes() - 2;
1040 int Qi =
m_base[1]->GetNumModes() - 2;
1041 int Ri =
m_base[2]->GetNumModes() - 2;
1045 return Pi * (2 * Qi - Pi - 1) / 2;
1049 return Pi * (2 * Ri - Pi - 1) / 2;
1053 return Qi * (2 * Ri - Qi - 1) / 2;
1059 ASSERTL2(i >= 0 && i <= 3,
"face id is out of range");
1063 return m_base[0]->GetNumPoints() *
m_base[1]->GetNumPoints();
1067 return m_base[0]->GetNumPoints() *
m_base[2]->GetNumPoints();
1071 return m_base[1]->GetNumPoints() *
m_base[2]->GetNumPoints();
1077 ASSERTL2((i >= 0) && (i <= 5),
"edge id is out of range");
1078 int P =
m_base[0]->GetNumModes();
1079 int Q =
m_base[1]->GetNumModes();
1080 int R =
m_base[2]->GetNumModes();
1086 else if (i == 1 || i == 2)
1099 ASSERTL2(i >= 0 && i <= 3,
"face id is out of range");
1100 ASSERTL2(j == 0 || j == 1,
"face direction is out of range");
1104 return m_base[j]->GetPointsKey();
1108 return m_base[2 * j]->GetPointsKey();
1112 return m_base[j + 1]->GetPointsKey();
1117 const std::vector<unsigned int> &nummodes,
int &modes_offset)
1120 nummodes[modes_offset], nummodes[modes_offset + 1],
1121 nummodes[modes_offset + 2]);
1131 ASSERTL2(i >= 0 && i <= 4,
"face id is out of range");
1132 ASSERTL2(k == 0 || k == 1,
"face direction out of range");
1165 for (
int k = 0; k < Qz; ++k)
1167 for (
int j = 0; j < Qy; ++j)
1169 for (
int i = 0; i < Qx; ++i)
1171 int s = i + Qx * (j + Qy * k);
1173 (eta_x[i] + 1.0) * (1.0 - eta_y[j]) * (1.0 - eta_z[k]) / 4 -
1175 xi_y[s] = (eta_y[j] + 1.0) * (1.0 - eta_z[k]) / 2 - 1.0;
1197 "Mapping not defined for this type of basis");
1200 if (useCoeffPacking ==
true)
1202 switch (localVertexId)
1226 ASSERTL0(
false,
"Vertex ID must be between 0 and 3");
1233 switch (localVertexId)
1257 ASSERTL0(
false,
"Vertex ID must be between 0 and 3");
1277 "BasisType is not a boundary interior form");
1280 "BasisType is not a boundary interior form");
1283 "BasisType is not a boundary interior form");
1285 int P =
m_base[0]->GetNumModes();
1286 int Q =
m_base[1]->GetNumModes();
1287 int R =
m_base[2]->GetNumModes();
1291 if (outarray.size() != nIntCoeffs)
1297 for (
int i = 2; i <
P; ++i)
1299 for (
int j = 1; j < Q - i - 1; ++j)
1301 for (
int k = 1; k < R - i - j; ++k)
1303 outarray[idx++] =
GetMode(i, j, k);
1316 "BasisType is not a boundary interior form");
1319 "BasisType is not a boundary interior form");
1322 "BasisType is not a boundary interior form");
1324 int P =
m_base[0]->GetNumModes();
1325 int Q =
m_base[1]->GetNumModes();
1326 int R =
m_base[2]->GetNumModes();
1333 if (outarray.size() != nBnd)
1338 for (i = 0; i <
P; ++i)
1343 for (j = 0; j < Q - i; j++)
1345 for (k = 0; k < R - i - j; ++k)
1347 outarray[idx++] =
GetMode(i, j, k);
1355 for (k = 0; k < R - i; ++k)
1357 outarray[idx++] =
GetMode(i, 0, k);
1359 for (j = 1; j < Q - i; ++j)
1361 outarray[idx++] =
GetMode(i, j, 0);
1371 int P = 0, Q = 0, idx = 0;
1372 int nFaceCoeffs = 0;
1378 Q =
m_base[1]->GetNumModes();
1382 Q =
m_base[2]->GetNumModes();
1387 Q =
m_base[2]->GetNumModes();
1390 ASSERTL0(
false,
"fid must be between 0 and 3");
1393 nFaceCoeffs =
P * (2 * Q -
P + 1) / 2;
1395 if (maparray.size() != nFaceCoeffs)
1404 for (i = 0; i <
P; ++i)
1406 for (j = 0; j < Q - i; ++j)
1408 maparray[idx++] =
GetMode(i, j, 0);
1414 for (i = 0; i <
P; ++i)
1416 for (k = 0; k < Q - i; ++k)
1418 maparray[idx++] =
GetMode(i, 0, k);
1424 for (j = 0; j <
P - 1; ++j)
1426 for (k = 0; k < Q - 1 - j; ++k)
1428 maparray[idx++] =
GetMode(1, j, k);
1430 if (j == 0 && k == 0)
1432 maparray[idx++] =
GetMode(0, 0, 1);
1434 if (j == 0 && k == Q - 2)
1436 for (
int r = 0; r < Q - 1; ++r)
1438 maparray[idx++] =
GetMode(0, 1, r);
1446 for (j = 0; j <
P; ++j)
1448 for (k = 0; k < Q - j; ++k)
1450 maparray[idx++] =
GetMode(0, j, k);
1455 ASSERTL0(
false,
"Element map not available.");
1464 int nummodesA = 0, nummodesB = 0, i, j, k, idx;
1467 "Method only implemented for Modified_A BasisType (x "
1468 "direction), Modified_B BasisType (y direction), and "
1469 "Modified_C BasisType(z direction)");
1471 int nFaceCoeffs = 0;
1476 nummodesA =
m_base[0]->GetNumModes();
1477 nummodesB =
m_base[1]->GetNumModes();
1480 nummodesA =
m_base[0]->GetNumModes();
1481 nummodesB =
m_base[2]->GetNumModes();
1485 nummodesA =
m_base[1]->GetNumModes();
1486 nummodesB =
m_base[2]->GetNumModes();
1489 ASSERTL0(
false,
"fid must be between 0 and 3");
1498 nFaceCoeffs =
P * (2 * Q -
P + 1) / 2;
1501 if (maparray.size() != nFaceCoeffs)
1506 if (signarray.size() != nFaceCoeffs)
1512 fill(signarray.data(), signarray.data() + nFaceCoeffs, 1);
1519 int minPA =
min(nummodesA,
P);
1520 int minQB =
min(nummodesB, Q);
1522 for (j = 0; j < minPA; ++j)
1525 for (k = 0; k < minQB - j; ++k, ++cnt)
1527 maparray[idx++] = cnt;
1530 cnt += nummodesB - minQB;
1532 for (k = nummodesB - j; k < Q - j; ++k)
1534 signarray[idx] = 0.0;
1535 maparray[idx++] = maparray[0];
1539 for (j = nummodesA; j <
P; ++j)
1541 for (k = 0; k < Q - j; ++k)
1543 signarray[idx] = 0.0;
1544 maparray[idx++] = maparray[0];
1551 for (i = 0; i <
P; ++i)
1553 for (j = 0; j < Q - i; ++j, idx++)
1557 signarray[idx] = (i % 2 ? -1 : 1);
1562 swap(maparray[0], maparray[Q]);
1564 for (i = 1; i < Q - 1; ++i)
1566 swap(maparray[i + 1], maparray[Q + i]);
1579 const int P =
m_base[0]->GetNumModes();
1580 const int Q =
m_base[1]->GetNumModes();
1581 const int R =
m_base[2]->GetNumModes();
1585 if (maparray.size() != nEdgeIntCoeffs)
1591 fill(maparray.data(), maparray.data() + nEdgeIntCoeffs, 0);
1594 if (signarray.size() != nEdgeIntCoeffs)
1600 fill(signarray.data(), signarray.data() + nEdgeIntCoeffs, 1);
1606 for (i = 0; i <
P - 2; ++i)
1608 maparray[i] =
GetMode(i + 2, 0, 0);
1612 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1619 for (i = 0; i < Q - 2; ++i)
1621 maparray[i] =
GetMode(1, i + 1, 0);
1625 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1632 for (i = 0; i < Q - 2; ++i)
1634 maparray[i] =
GetMode(0, i + 2, 0);
1638 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1645 for (i = 0; i < R - 2; ++i)
1647 maparray[i] =
GetMode(0, 0, i + 2);
1651 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1658 for (i = 0; i < R - 2; ++i)
1660 maparray[i] =
GetMode(1, 0, i + 1);
1664 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1671 for (i = 0; i < R - 2; ++i)
1673 maparray[i] =
GetMode(0, 1, i + 1);
1677 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1684 ASSERTL0(
false,
"Edge not defined.");
1694 const int P =
m_base[0]->GetNumModes();
1695 const int Q =
m_base[1]->GetNumModes();
1696 const int R =
m_base[2]->GetNumModes();
1700 if (maparray.size() != nFaceIntCoeffs)
1705 if (signarray.size() != nFaceIntCoeffs)
1711 fill(signarray.data(), signarray.data() + nFaceIntCoeffs, 1);
1718 for (i = 2; i <
P; ++i)
1720 for (j = 1; j < Q - i; ++j)
1722 if ((
int)faceOrient == 7)
1724 signarray[idx] = (i % 2 ? -1 : 1);
1726 maparray[idx++] =
GetMode(i, j, 0);
1732 for (i = 2; i <
P; ++i)
1734 for (k = 1; k < R - i; ++k)
1736 if ((
int)faceOrient == 7)
1738 signarray[idx] = (i % 2 ? -1 : 1);
1740 maparray[idx++] =
GetMode(i, 0, k);
1746 for (j = 1; j < Q - 1; ++j)
1748 for (k = 1; k < R - 1 - j; ++k)
1750 if ((
int)faceOrient == 7)
1752 signarray[idx] = ((j + 1) % 2 ? -1 : 1);
1754 maparray[idx++] =
GetMode(1, j, k);
1760 for (j = 2; j < Q; ++j)
1762 for (k = 1; k < R - j; ++k)
1764 if ((
int)faceOrient == 7)
1766 signarray[idx] = (j % 2 ? -1 : 1);
1768 maparray[idx++] =
GetMode(0, j, k);
1773 ASSERTL0(
false,
"Face interior map not available.");
1791 int nq0 =
m_base[0]->GetNumPoints();
1792 int nq1 =
m_base[1]->GetNumPoints();
1793 int nq2 =
m_base[2]->GetNumPoints();
1803 nq =
max(nq0,
max(nq1, nq2));
1817 for (
int i = 0; i < nq; ++i)
1819 for (
int j = 0; j < nq - i; ++j)
1821 for (
int k = 0; k < nq - i - j; ++k, ++cnt)
1824 coords[cnt][0] = -1.0 + 2 * k / (
NekDouble)(nq - 1);
1825 coords[cnt][1] = -1.0 + 2 * j / (
NekDouble)(nq - 1);
1826 coords[cnt][2] = -1.0 + 2 * i / (
NekDouble)(nq - 1);
1831 for (
int i = 0; i < neq; ++i)
1835 I[0] =
m_base[0]->GetI(coll);
1836 I[1] =
m_base[1]->GetI(coll + 1);
1837 I[2] =
m_base[2]->GetI(coll + 2);
1841 for (
int k = 0; k < nq2; ++k)
1843 for (
int j = 0; j < nq1; ++j)
1846 fac = (I[1]->GetPtr())[j] * (I[2]->GetPtr())[k];
1850 Mat->GetRawPtr() + k * nq0 * nq1 * neq +
1914 const int Q =
m_base[1]->GetNumModes();
1915 const int R =
m_base[2]->GetNumModes();
1917 int i, j, q_hat, k_hat;
1921 for (i = 0; i < I; ++i)
1927 cnt += q_hat * (q_hat + 1) / 2 + k_hat * (Q - i);
1932 for (j = 0; j < J; ++j)
1950 int nquad0 =
m_base[0]->GetNumPoints();
1951 int nquad1 =
m_base[1]->GetNumPoints();
1952 int nquad2 =
m_base[2]->GetNumPoints();
1962 for (i = 0; i < nquad1 * nquad2; ++i)
1965 1, &outarray[0] + i * nquad0, 1);
1971 case LibUtilities::eGaussRadauMAlpha1Beta0:
1972 for (j = 0; j < nquad2; ++j)
1974 for (i = 0; i < nquad1; ++i)
1977 &outarray[0] + i * nquad0 + j * nquad0 * nquad1,
1984 for (j = 0; j < nquad2; ++j)
1986 for (i = 0; i < nquad1; ++i)
1989 &outarray[0] + i * nquad0 + j * nquad0 * nquad1,
1999 case LibUtilities::eGaussRadauMAlpha2Beta0:
2000 for (i = 0; i < nquad2; ++i)
2003 &outarray[0] + i * nquad0 * nquad1, 1);
2007 case LibUtilities::eGaussRadauMAlpha1Beta0:
2008 for (i = 0; i < nquad2; ++i)
2010 Blas::Dscal(nquad0 * nquad1, 0.25 * (1 - z2[i]) * w2[i],
2011 &outarray[0] + i * nquad0 * nquad1, 1);
2015 for (i = 0; i < nquad2; ++i)
2018 0.25 * (1 - z2[i]) * (1 - z2[i]) * w2[i],
2019 &outarray[0] + i * nquad0 * nquad1, 1);
2036 int qa =
m_base[0]->GetNumPoints();
2037 int qb =
m_base[1]->GetNumPoints();
2038 int qc =
m_base[2]->GetNumPoints();
2039 int nmodes_a =
m_base[0]->GetNumModes();
2040 int nmodes_b =
m_base[1]->GetNumModes();
2041 int nmodes_c =
m_base[2]->GetNumModes();
2055 int i, j, k, cnt = 0;
2058 OrthoExp.
FwdTrans(array, orthocoeffs);
2068 for (i = 0; i < nmodes_a; ++i)
2070 for (j = 0; j < nmodes_b - j; ++j)
2073 pow((1.0 * i) / (nmodes_a - 1), cutoff * nmodes_a),
2074 pow((1.0 * j) / (nmodes_b - 1), cutoff * nmodes_b));
2076 for (k = 0; k < nmodes_c - i - j; ++k)
2079 std::max(fac1, pow((1.0 * k) / (nmodes_c - 1),
2080 cutoff * nmodes_c));
2082 orthocoeffs[cnt] *= SvvDiffCoeff * fac;
2098 max_abc =
max(max_abc, 0);
2101 for (i = 0; i < nmodes_a; ++i)
2103 for (j = 0; j < nmodes_b - j; ++j)
2105 int maxij =
max(i, j);
2107 for (k = 0; k < nmodes_c - i - j; ++k)
2109 int maxijk =
max(maxij, k);
2132 int cutoff_a = (int)(SVVCutOff * nmodes_a);
2133 int cutoff_b = (int)(SVVCutOff * nmodes_b);
2134 int cutoff_c = (int)(SVVCutOff * nmodes_c);
2135 int nmodes =
min(
min(nmodes_a, nmodes_b), nmodes_c);
2140 for (i = 0; i < nmodes_a; ++i)
2142 for (j = 0; j < nmodes_b - i; ++j)
2144 for (k = 0; k < nmodes_c - i - j; ++k)
2146 if (i + j + k >= cutoff)
2148 orthocoeffs[cnt] *= ((SvvDiffCoeff)*exp(
2149 -(i + j + k - nmodes) * (i + j + k - nmodes) /
2150 ((
NekDouble)((i + j + k - cutoff + epsilon) *
2151 (i + j + k - cutoff + epsilon)))));
2155 orthocoeffs[cnt] *= 0.0;
2164 OrthoExp.
BwdTrans(orthocoeffs, array);
2171 int nquad0 =
m_base[0]->GetNumPoints();
2172 int nquad1 =
m_base[1]->GetNumPoints();
2173 int nquad2 =
m_base[2]->GetNumPoints();
2174 int nqtot = nquad0 * nquad1 * nquad2;
2175 int nmodes0 =
m_base[0]->GetNumModes();
2176 int nmodes1 =
m_base[1]->GetNumModes();
2177 int nmodes2 =
m_base[2]->GetNumModes();
2178 int numMax = nmodes0;
2200 bortho0, bortho1, bortho2);
2203 OrthoTetExp->FwdTrans(phys_tmp, coeff);
2209 for (
int u = 0; u < numMin; ++u)
2211 for (
int i = 0; i < numMin - u; ++i)
2214 tmp2 = coeff_tmp1 + cnt, 1);
2215 cnt += numMax - u - i;
2217 for (
int i = numMin; i < numMax - u; ++i)
2219 cnt += numMax - u - i;
2223 OrthoTetExp->BwdTrans(coeff_tmp1, phys_tmp);
2230 int np0 =
m_base[0]->GetNumPoints();
2231 int np1 =
m_base[1]->GetNumPoints();
2232 int np2 =
m_base[2]->GetNumPoints();
2233 int np =
max(np0,
max(np1, np2));
2244 for (
int i = 0; i < np - 1; ++i)
2246 planep1 += (np - i) * (np - i + 1) / 2;
2251 for (
int j = 0; j < np - i - 1; ++j)
2253 rowp1 += np - i - j;
2254 row1p1 += np - i - j - 1;
2255 for (
int k = 0; k < np - i - j - 2; ++k)
2257 conn[cnt++] = plane + row + k + 1;
2258 conn[cnt++] = plane + row + k;
2259 conn[cnt++] = plane + rowp1 + k;
2260 conn[cnt++] = planep1 + row1 + k;
2262 conn[cnt++] = plane + row + k + 1;
2263 conn[cnt++] = plane + rowp1 + k + 1;
2264 conn[cnt++] = planep1 + row1 + k + 1;
2265 conn[cnt++] = planep1 + row1 + k;
2267 conn[cnt++] = plane + rowp1 + k + 1;
2268 conn[cnt++] = plane + row + k + 1;
2269 conn[cnt++] = plane + rowp1 + k;
2270 conn[cnt++] = planep1 + row1 + k;
2272 conn[cnt++] = planep1 + row1 + k;
2273 conn[cnt++] = planep1 + row1p1 + k;
2274 conn[cnt++] = plane + rowp1 + k;
2275 conn[cnt++] = plane + rowp1 + k + 1;
2277 conn[cnt++] = planep1 + row1 + k;
2278 conn[cnt++] = planep1 + row1p1 + k;
2279 conn[cnt++] = planep1 + row1 + k + 1;
2280 conn[cnt++] = plane + rowp1 + k + 1;
2282 if (k < np - i - j - 3)
2284 conn[cnt++] = plane + rowp1 + k + 1;
2285 conn[cnt++] = planep1 + row1p1 + k + 1;
2286 conn[cnt++] = planep1 + row1 + k + 1;
2287 conn[cnt++] = planep1 + row1p1 + k;
2291 conn[cnt++] = plane + row + np - i - j - 1;
2292 conn[cnt++] = plane + row + np - i - j - 2;
2293 conn[cnt++] = plane + rowp1 + np - i - j - 2;
2294 conn[cnt++] = planep1 + row1 + np - i - j - 2;
2297 row1 += np - i - j - 1;
2299 plane += (np - i) * (np - i + 1) / 2;
#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.
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Backward transformation from coefficient space to physical space.
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Forward transformation from physical space to coefficient space.
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
void PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Array< OneD, LibUtilities::BasisSharedPtr > m_base
MatrixType GetMatrixType() const
NekDouble GetConstFactor(const ConstFactorType &factor) const
bool ConstFactorExists(const ConstFactorType &factor) const
int v_GetNtraces() const override
int v_GetTraceNcoeffs(const int i) const override
NekDouble v_PhysEvalFirstDeriv(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
void v_GetTraceNumModes(const int fid, int &numModes0, int &numModes1, Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2) override
int v_GetTraceIntNcoeffs(const int i) const override
void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray) override
void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_GetInteriorMap(Array< OneD, unsigned int > &outarray) override
int v_GetNedges() const override
int GetMode(const int i, const int j, const int k)
Compute the mode number in the expansion for a particular tensorial combination.
int v_GetEdgeNcoeffs(const int i) const override
bool v_IsBoundaryInteriorExpansion() const override
LibUtilities::PointsKey v_GetTracePointsKey(const int i, const int j) const override
NekDouble v_PhysEvaluateBasis(const Array< OneD, const NekDouble > &coords, int mode) final
void v_GetEdgeInteriorToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2) override
void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
StdTetExp(const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
int v_GetTraceNumPoints(const int i) const override
void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_BwdTrans_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2) override
void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey) override
DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey) override
void v_LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta) override
void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray) override
void v_GetSimplexEquiSpacedConnectivity(Array< OneD, int > &conn, bool standard=true) override
void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
int v_NumBndryCoeffs() const override
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_StdPhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2) override
void v_GetCoords(Array< OneD, NekDouble > &coords_x, Array< OneD, NekDouble > &coords_y, Array< OneD, NekDouble > &coords_z) override
void v_GetTraceCoeffMap(const unsigned int fid, Array< OneD, unsigned int > &maparray) override
void v_GetTraceInteriorToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2) override
void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_GetElmtTraceToTraceMap(const unsigned int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1) override
int v_GetVertexMap(int localVertexId, bool useCoeffPacking=false) override
int v_GetNverts() const override
void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dx, Array< OneD, NekDouble > &out_dy, Array< OneD, NekDouble > &out_dz) override
Calculate the derivative of the physical points.
const LibUtilities::BasisKey v_GetTraceBasisKey(const int i, const int k, bool UseGLL=false) const override
LibUtilities::ShapeType v_DetShapeType() const override
void v_LocCollapsedToLocCoord(const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi) override
LibUtilities::ShapeType DetShapeType() const
int v_NumDGBndryCoeffs() const override
void v_IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2) override
int v_CalcNumberOfCoefficients(const std::vector< unsigned int > &nummodes, int &modes_offset) override
DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override
static void Dgemv(const char &trans, const int &m, const int &n, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = alpha A x plus beta y where A[m x n].
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 getNumberOfBndCoefficients(int Na, int Nb, int Nc)
int getNumberOfCoefficients(int Na, int Nb, int Nc)
@ 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
@ eFactorSVVDGKerDiffCoeff
@ eFactorSVVPowerKerDiffCoeff
std::shared_ptr< StdTetExp > StdTetExpSharedPtr
LibUtilities::BasisKey EvaluateTriFaceBasisKey(const int facedir, const LibUtilities::BasisSharedPtr &faceDirBasis, bool UseGLL)
const int kSVVDGFiltermodesmin
const int kSVVDGFiltermodesmax
const NekDouble kSVVDGFilter[9][11]
@ ePhysInterpToEquiSpaced
@ eDir1BwdDir1_Dir2FwdDir2
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 Zero(int n, T *x, const int incx)
Zero vector.
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector 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)