63 : StdExpansion(Ba.GetNumModes() * Bb.GetNumModes() * Bc.GetNumModes(), 3,
65 StdExpansion3D(Ba.GetNumModes() * Bb.GetNumModes() * Bc.GetNumModes(), Ba,
69 std::bind(&
Expansion3D::CreateMatrix, this, std::placeholders::_1),
70 std::string(
"HexExpMatrix")),
71 m_staticCondMatrixManager(std::bind(&
Expansion::CreateStaticCondMatrix,
72 this, std::placeholders::_1),
73 std::string(
"HexExpStaticCondMatrix"))
83 : StdExpansion(T), StdExpansion3D(T), StdHexExp(T),
Expansion(T),
85 m_staticCondMatrixManager(T.m_staticCondMatrixManager)
105 int nquad0 =
m_base[0]->GetNumPoints();
106 int nquad1 =
m_base[1]->GetNumPoints();
107 int nquad2 =
m_base[2]->GetNumPoints();
117 (
NekDouble *)&inarray[0], 1, &tmp[0], 1);
122 (
NekDouble *)&inarray[0], 1, &tmp[0], 1);
126 returnVal = StdHexExp::v_Integral(tmp);
149 int nquad0 =
m_base[0]->GetNumPoints();
150 int nquad1 =
m_base[1]->GetNumPoints();
151 int nquad2 =
m_base[2]->GetNumPoints();
152 int ntot = nquad0 * nquad1 * nquad2;
160 StdHexExp::v_PhysDeriv(inarray, Diff0, Diff1, Diff2);
166 Vmath::Vmul(ntot, &df[0][0], 1, &Diff0[0], 1, &out_d0[0], 1);
167 Vmath::Vvtvp(ntot, &df[1][0], 1, &Diff1[0], 1, &out_d0[0], 1,
169 Vmath::Vvtvp(ntot, &df[2][0], 1, &Diff2[0], 1, &out_d0[0], 1,
175 Vmath::Vmul(ntot, &df[3][0], 1, &Diff0[0], 1, &out_d1[0], 1);
176 Vmath::Vvtvp(ntot, &df[4][0], 1, &Diff1[0], 1, &out_d1[0], 1,
178 Vmath::Vvtvp(ntot, &df[5][0], 1, &Diff2[0], 1, &out_d1[0], 1,
184 Vmath::Vmul(ntot, &df[6][0], 1, &Diff0[0], 1, &out_d2[0], 1);
185 Vmath::Vvtvp(ntot, &df[7][0], 1, &Diff1[0], 1, &out_d2[0], 1,
187 Vmath::Vvtvp(ntot, &df[8][0], 1, &Diff2[0], 1, &out_d2[0], 1,
195 Vmath::Smul(ntot, df[0][0], &Diff0[0], 1, &out_d0[0], 1);
196 Blas::Daxpy(ntot, df[1][0], &Diff1[0], 1, &out_d0[0], 1);
197 Blas::Daxpy(ntot, df[2][0], &Diff2[0], 1, &out_d0[0], 1);
202 Vmath::Smul(ntot, df[3][0], &Diff0[0], 1, &out_d1[0], 1);
203 Blas::Daxpy(ntot, df[4][0], &Diff1[0], 1, &out_d1[0], 1);
204 Blas::Daxpy(ntot, df[5][0], &Diff2[0], 1, &out_d1[0], 1);
209 Vmath::Smul(ntot, df[6][0], &Diff0[0], 1, &out_d2[0], 1);
210 Blas::Daxpy(ntot, df[7][0], &Diff1[0], 1, &out_d2[0], 1);
211 Blas::Daxpy(ntot, df[8][0], &Diff2[0], 1, &out_d2[0], 1);
251 ASSERTL1(
false,
"input dir is out of range");
264 int nquad0 =
m_base[0]->GetNumPoints();
265 int nquad1 =
m_base[1]->GetNumPoints();
266 int nquad2 =
m_base[2]->GetNumPoints();
267 int ntot = nquad0 * nquad1 * nquad2;
275 StdHexExp::v_PhysDeriv(inarray, Diff0, Diff1, Diff2);
280 Vmath::Vmul(ntot, &dfdir[0][0], 1, &Diff0[0], 1, &outarray[0], 1);
281 Vmath::Vvtvp(ntot, &dfdir[1][0], 1, &Diff1[0], 1, &outarray[0], 1,
283 Vmath::Vvtvp(ntot, &dfdir[2][0], 1, &Diff2[0], 1, &outarray[0], 1,
302 if (
m_base[0]->Collocation() &&
m_base[1]->Collocation() &&
319 out = (*matsys) * in;
376 int nquad0 =
m_base[0]->GetNumPoints();
377 int nquad1 =
m_base[1]->GetNumPoints();
378 int nquad2 =
m_base[2]->GetNumPoints();
379 int order0 =
m_base[0]->GetNumModes();
380 int order1 =
m_base[1]->GetNumModes();
383 order0 * order1 * nquad2);
385 if (multiplybyweights)
392 tmp, outarray, wsp,
true,
true,
true);
398 inarray, outarray, wsp,
true,
true,
true);
433 ASSERTL1((dir == 0) || (dir == 1) || (dir == 2),
"Invalid direction.");
435 const int nq0 =
m_base[0]->GetNumPoints();
436 const int nq1 =
m_base[1]->GetNumPoints();
437 const int nq2 =
m_base[2]->GetNumPoints();
438 const int nq = nq0 * nq1 * nq2;
439 const int nm0 =
m_base[0]->GetNumModes();
440 const int nm1 =
m_base[1]->GetNumModes();
460 m_base[2]->GetBdata(), tmp2, outarray, wsp,
464 m_base[2]->GetBdata(), tmp3, tmp5, wsp,
true,
469 m_base[2]->GetDbdata(), tmp4, tmp5, wsp,
true,
478 ASSERTL1((dir == 0) || (dir == 1) || (dir == 2),
"Invalid direction.");
480 const int nq0 =
m_base[0]->GetNumPoints();
481 const int nq1 =
m_base[1]->GetNumPoints();
482 const int nq2 =
m_base[2]->GetNumPoints();
483 const int nq = nq0 * nq1 * nq2;
498 Vmath::Vmul(nq, &df[3 * dir][0], 1, tmp1.get(), 1, tmp2.get(), 1);
499 Vmath::Vmul(nq, &df[3 * dir + 1][0], 1, tmp1.get(), 1, tmp3.get(), 1);
500 Vmath::Vmul(nq, &df[3 * dir + 2][0], 1, tmp1.get(), 1, tmp4.get(), 1);
504 Vmath::Smul(nq, df[3 * dir][0], tmp1.get(), 1, tmp2.get(), 1);
505 Vmath::Smul(nq, df[3 * dir + 1][0], tmp1.get(), 1, tmp3.get(), 1);
506 Vmath::Smul(nq, df[3 * dir + 2][0], tmp1.get(), 1, tmp4.get(), 1);
522 const int nq0 =
m_base[0]->GetNumPoints();
523 const int nq1 =
m_base[1]->GetNumPoints();
524 const int nq2 =
m_base[2]->GetNumPoints();
525 const int nq = nq0 * nq1 * nq2;
526 const int nm0 =
m_base[0]->GetNumModes();
527 const int nm1 =
m_base[1]->GetNumModes();
545 Vmath::Vmul(nq, &dfdir[0][0], 1, tmp1.get(), 1, tmp2.get(), 1);
546 Vmath::Vmul(nq, &dfdir[1][0], 1, tmp1.get(), 1, tmp3.get(), 1);
547 Vmath::Vmul(nq, &dfdir[2][0], 1, tmp1.get(), 1, tmp4.get(), 1);
550 m_base[2]->GetBdata(), tmp2, outarray, wsp,
554 m_base[2]->GetBdata(), tmp3, tmp5, wsp,
true,
560 m_base[2]->GetDbdata(), tmp4, tmp5, wsp,
true,
580 return StdExpansion3D::v_PhysEvaluate(Lcoord, physvals);
589 m_geom->GetLocCoords(coord, Lcoord);
590 return StdExpansion3D::v_PhysEvaluate(Lcoord, physvals);
595 std::array<NekDouble, 3> &firstOrderDerivs)
599 m_geom->GetLocCoords(coord, Lcoord);
600 return StdHexExp::v_PhysEvaluate(Lcoord, inarray, firstOrderDerivs);
607 m_base[2]->GetBasisKey());
613 m_base[0]->GetPointsKey());
615 m_base[1]->GetPointsKey());
617 m_base[2]->GetPointsKey());
635 ASSERTL1(Lcoords[0] >= -1.0 && Lcoords[0] <= 1.0 && Lcoords[1] >= -1.0 &&
636 Lcoords[1] <= 1.0 && Lcoords[2] >= -1.0 && Lcoords[2] <= 1.0,
637 "Local coordinates are not in region [-1,1]");
641 for (i = 0; i <
m_geom->GetCoordim(); ++i)
643 coords[i] =
m_geom->GetCoord(i, Lcoords);
665 const NekDouble *data,
const std::vector<unsigned int> &nummodes,
666 const int mode_offset,
NekDouble *coeffs,
667 std::vector<LibUtilities::BasisType> &fromType)
669 int data_order0 = nummodes[mode_offset];
670 int fillorder0 = min(
m_base[0]->GetNumModes(), data_order0);
671 int data_order1 = nummodes[mode_offset + 1];
672 int order1 =
m_base[1]->GetNumModes();
673 int fillorder1 = min(order1, data_order1);
674 int data_order2 = nummodes[mode_offset + 2];
675 int order2 =
m_base[2]->GetNumModes();
676 int fillorder2 = min(order2, data_order2);
688 m_base[0]->GetPointsKey()),
690 m_base[1]->GetPointsKey()),
692 m_base[2]->GetPointsKey()));
695 m_base[2]->GetBasisKey());
717 "Extraction routine not set up for this basis");
719 "Extraction routine not set up for this basis");
722 for (j = 0; j < fillorder0; ++j)
724 for (i = 0; i < fillorder1; ++i)
726 Vmath::Vcopy(fillorder2, &data[cnt], 1, &coeffs[cnt1], 1);
732 for (i = fillorder1; i < data_order1; ++i)
737 for (i = fillorder1; i < order1; ++i)
762 ASSERTL0(
false,
"basis is either not set up or not "
769 int nquad0 =
m_base[0]->GetNumPoints();
770 int nquad1 =
m_base[1]->GetNumPoints();
771 int nquad2 =
m_base[2]->GetNumPoints();
783 if (outarray.size() != nq0 * nq1)
788 for (
int i = 0; i < nquad0 * nquad1; ++i)
798 if (outarray.size() != nq0 * nq1)
804 for (
int k = 0; k < nquad2; k++)
806 for (
int i = 0; i < nquad0; ++i)
808 outarray[k * nquad0 + i] = nquad0 * nquad1 * k + i;
817 if (outarray.size() != nq0 * nq1)
822 for (
int i = 0; i < nquad1 * nquad2; i++)
824 outarray[i] = nquad0 - 1 + i * nquad0;
832 if (outarray.size() != nq0 * nq1)
837 for (
int k = 0; k < nquad2; k++)
839 for (
int i = 0; i < nquad0; i++)
841 outarray[k * nquad0 + i] =
842 (nquad0 * (nquad1 - 1)) + (k * nquad0 * nquad1) + i;
851 if (outarray.size() != nq0 * nq1)
856 for (
int i = 0; i < nquad1 * nquad2; i++)
858 outarray[i] = i * nquad0;
865 if (outarray.size() != nq0 * nq1)
870 for (
int i = 0; i < nquad0 * nquad1; i++)
872 outarray[i] = nquad0 * nquad1 * (nquad2 - 1) + i;
877 ASSERTL0(
false,
"face value (> 5) is out of range");
890 for (i = 0; i < ptsKeys.size(); ++i)
901 geomFactors->GetDerivFactors(ptsKeys);
914 for (i = 0; i < vCoordDim; ++i)
919 size_t nqb = nq_face;
933 for (i = 0; i < vCoordDim; ++i)
935 normal[i][0] = -df[3 * i + 2][0];
939 for (i = 0; i < vCoordDim; ++i)
941 normal[i][0] = -df[3 * i + 1][0];
945 for (i = 0; i < vCoordDim; ++i)
947 normal[i][0] = df[3 * i][0];
951 for (i = 0; i < vCoordDim; ++i)
953 normal[i][0] = df[3 * i + 1][0];
957 for (i = 0; i < vCoordDim; ++i)
959 normal[i][0] = -df[3 * i][0];
963 for (i = 0; i < vCoordDim; ++i)
965 normal[i][0] = df[3 * i + 2][0];
969 ASSERTL0(
false,
"face is out of range (edge < 5)");
974 for (i = 0; i < vCoordDim; ++i)
976 fac += normal[i][0] * normal[i][0];
978 fac = 1.0 /
sqrt(fac);
981 for (i = 0; i < vCoordDim; ++i)
983 Vmath::Fill(nq_face, fac * normal[i][0], normal[i], 1);
990 int nqe0 = ptsKeys[0].GetNumPoints();
991 int nqe1 = ptsKeys[1].GetNumPoints();
992 int nqe2 = ptsKeys[2].GetNumPoints();
993 int nqe01 = nqe0 * nqe1;
994 int nqe02 = nqe0 * nqe2;
995 int nqe12 = nqe1 * nqe2;
998 if (face == 0 || face == 5)
1002 else if (face == 1 || face == 3)
1023 for (j = 0; j < nqe; ++j)
1025 normals[j] = -df[2][j] * jac[j];
1026 normals[nqe + j] = -df[5][j] * jac[j];
1027 normals[2 * nqe + j] = -df[8][j] * jac[j];
1028 faceJac[j] = jac[j];
1031 points0 = ptsKeys[0];
1032 points1 = ptsKeys[1];
1035 for (j = 0; j < nqe0; ++j)
1037 for (k = 0; k < nqe2; ++k)
1039 int idx = j + nqe01 * k;
1040 normals[j + k * nqe0] = -df[1][idx] * jac[idx];
1041 normals[nqe + j + k * nqe0] = -df[4][idx] * jac[idx];
1042 normals[2 * nqe + j + k * nqe0] =
1043 -df[7][idx] * jac[idx];
1044 faceJac[j + k * nqe0] = jac[idx];
1047 points0 = ptsKeys[0];
1048 points1 = ptsKeys[2];
1051 for (j = 0; j < nqe1; ++j)
1053 for (k = 0; k < nqe2; ++k)
1055 int idx = nqe0 - 1 + nqe0 * j + nqe01 * k;
1056 normals[j + k * nqe1] = df[0][idx] * jac[idx];
1057 normals[nqe + j + k * nqe1] = df[3][idx] * jac[idx];
1058 normals[2 * nqe + j + k * nqe1] = df[6][idx] * jac[idx];
1059 faceJac[j + k * nqe1] = jac[idx];
1062 points0 = ptsKeys[1];
1063 points1 = ptsKeys[2];
1066 for (j = 0; j < nqe0; ++j)
1068 for (k = 0; k < nqe2; ++k)
1070 int idx = nqe0 * (nqe1 - 1) + j + nqe01 * k;
1071 normals[j + k * nqe0] = df[1][idx] * jac[idx];
1072 normals[nqe + j + k * nqe0] = df[4][idx] * jac[idx];
1073 normals[2 * nqe + j + k * nqe0] = df[7][idx] * jac[idx];
1074 faceJac[j + k * nqe0] = jac[idx];
1077 points0 = ptsKeys[0];
1078 points1 = ptsKeys[2];
1081 for (j = 0; j < nqe1; ++j)
1083 for (k = 0; k < nqe2; ++k)
1085 int idx = j * nqe0 + nqe01 * k;
1086 normals[j + k * nqe1] = -df[0][idx] * jac[idx];
1087 normals[nqe + j + k * nqe1] = -df[3][idx] * jac[idx];
1088 normals[2 * nqe + j + k * nqe1] =
1089 -df[6][idx] * jac[idx];
1090 faceJac[j + k * nqe1] = jac[idx];
1093 points0 = ptsKeys[1];
1094 points1 = ptsKeys[2];
1097 for (j = 0; j < nqe01; ++j)
1099 int idx = j + nqe01 * (nqe2 - 1);
1100 normals[j] = df[2][idx] * jac[idx];
1101 normals[nqe + j] = df[5][idx] * jac[idx];
1102 normals[2 * nqe + j] = df[8][idx] * jac[idx];
1103 faceJac[j] = jac[idx];
1105 points0 = ptsKeys[0];
1106 points1 = ptsKeys[1];
1109 ASSERTL0(
false,
"face is out of range (face < 5)");
1118 Vmath::Sdiv(nq_face, 1.0, &work[0], 1, &work[0], 1);
1126 Vmath::Vmul(nq_face, work, 1, normal[i], 1, normal[i], 1);
1133 Vmath::Vvtvp(nq_face, normal[i], 1, normal[i], 1, work, 1, work, 1);
1143 Vmath::Vmul(nq_face, normal[i], 1, work, 1, normal[i], 1);
1155 StdExpansion::MassMatrixOp_MatFree(inarray, outarray, mkey);
1170 StdExpansion::LaplacianMatrixOp_MatFree(k1, k2, inarray, outarray, mkey);
1178 StdExpansion::WeakDerivMatrixOp_MatFree(i, inarray, outarray, mkey);
1185 StdExpansion::WeakDirectionalDerivMatrixOp_MatFree(inarray, outarray, mkey);
1192 StdExpansion::MassLevelCurvatureMatrixOp_MatFree(inarray, outarray, mkey);
1217 int n_coeffs = inarray.size();
1218 int nmodes0 =
m_base[0]->GetNumModes();
1219 int nmodes1 =
m_base[1]->GetNumModes();
1220 int nmodes2 =
m_base[2]->GetNumModes();
1221 int numMax = nmodes0;
1249 int cnt = 0, cnt2 = 0;
1251 for (
int u = 0; u < numMin + 1; ++u)
1253 for (
int i = 0; i < numMin; ++i)
1256 tmp2 = coeff_tmp1 + cnt, 1);
1262 tmp4 = coeff_tmp2 + cnt2, 1);
1264 cnt2 = u * nmodes0 * nmodes1;
1292 StdHexExp::v_SVVLaplacianFilter(array, mkey);
1317 returnval = StdHexExp::v_GenMatrix(mkey);
1332 return tmp->GetStdMatrix(mkey);
1366 int nquad0 =
m_base[0]->GetNumPoints();
1367 int nquad1 =
m_base[1]->GetNumPoints();
1368 int nquad2 =
m_base[2]->GetNumPoints();
1369 int nqtot = nquad0 * nquad1 * nquad2;
1371 ASSERTL1(wsp.size() >= 6 * nqtot,
"Insufficient workspace size.");
1400 StdExpansion3D::PhysTensorDeriv(inarray, wsp0, wsp1, wsp2);
1406 Vmath::Vvtvvtp(nqtot, &metric00[0], 1, &wsp0[0], 1, &metric01[0], 1,
1407 &wsp1[0], 1, &wsp3[0], 1);
1408 Vmath::Vvtvp(nqtot, &metric02[0], 1, &wsp2[0], 1, &wsp3[0], 1, &wsp3[0], 1);
1409 Vmath::Vvtvvtp(nqtot, &metric01[0], 1, &wsp0[0], 1, &metric11[0], 1,
1410 &wsp1[0], 1, &wsp4[0], 1);
1411 Vmath::Vvtvp(nqtot, &metric12[0], 1, &wsp2[0], 1, &wsp4[0], 1, &wsp4[0], 1);
1412 Vmath::Vvtvvtp(nqtot, &metric02[0], 1, &wsp0[0], 1, &metric12[0], 1,
1413 &wsp1[0], 1, &wsp5[0], 1);
1414 Vmath::Vvtvp(nqtot, &metric22[0], 1, &wsp2[0], 1, &wsp5[0], 1, &wsp5[0], 1);
1437 const unsigned int dim = 3;
1443 for (
unsigned int i = 0; i < dim; ++i)
1445 for (
unsigned int j = i; j < dim; ++j)
1482 if (d0factors.size() != 6)
1489 if (d0factors[0].size() != nquad0 * nquad1)
1499 if (d0factors[1].size() != nquad0 * nquad2)
1509 if (d0factors[2].size() != nquad1 * nquad2)
1533 int ncoords = normal_0.size();
1538 for (
int i = 0; i < nquad0 * nquad1; ++i)
1540 d0factors[0][i] = df[0][i] * normal_0[0][i];
1541 d1factors[0][i] = df[1][i] * normal_0[0][i];
1542 d2factors[0][i] = df[2][i] * normal_0[0][i];
1545 df[0][nquad0 * nquad1 * (nquad2 - 1) + i] * normal_5[0][i];
1547 df[1][nquad0 * nquad1 * (nquad2 - 1) + i] * normal_5[0][i];
1549 df[2][nquad0 * nquad1 * (nquad2 - 1) + i] * normal_5[0][i];
1552 for (
int n = 1; n < ncoords; ++n)
1554 for (
int i = 0; i < nquad0 * nquad1; ++i)
1556 d0factors[0][i] += df[3 * n][i] * normal_0[n][i];
1557 d1factors[0][i] += df[3 * n + 1][i] * normal_0[n][i];
1558 d2factors[0][i] += df[3 * n + 2][i] * normal_0[n][i];
1561 df[3 * n][nquad0 * nquad1 * (nquad2 - 1) + i] *
1564 df[3 * n + 1][nquad0 * nquad1 * (nquad2 - 1) + i] *
1567 df[3 * n + 2][nquad0 * nquad1 * (nquad2 - 1) + i] *
1573 for (
int j = 0; j < nquad2; ++j)
1575 for (
int i = 0; i < nquad0; ++i)
1577 d0factors[1][j * nquad0 + i] = df[0][j * nquad0 * nquad1 + i] *
1578 normal_1[0][j * nquad0 + i];
1579 d1factors[1][j * nquad0 + i] = df[1][j * nquad0 * nquad1 + i] *
1580 normal_1[0][j * nquad0 + i];
1581 d2factors[1][j * nquad0 + i] = df[2][j * nquad0 * nquad1 + i] *
1582 normal_1[0][j * nquad0 + i];
1584 d0factors[3][j * nquad0 + i] =
1585 df[0][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1586 normal_3[0][j * nquad0 + i];
1587 d1factors[3][j * nquad0 + i] =
1588 df[1][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1589 normal_3[0][j * nquad0 + i];
1590 d2factors[3][j * nquad0 + i] =
1591 df[2][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1592 normal_3[0][j * nquad0 + i];
1596 for (
int n = 1; n < ncoords; ++n)
1598 for (
int j = 0; j < nquad2; ++j)
1600 for (
int i = 0; i < nquad0; ++i)
1602 d0factors[1][j * nquad0 + i] +=
1603 df[3 * n][j * nquad0 * nquad1 + i] *
1604 normal_1[0][j * nquad0 + i];
1605 d1factors[1][j * nquad0 + i] +=
1606 df[3 * n + 1][j * nquad0 * nquad1 + i] *
1607 normal_1[0][j * nquad0 + i];
1608 d2factors[1][j * nquad0 + i] +=
1609 df[3 * n + 2][j * nquad0 * nquad1 + i] *
1610 normal_1[0][j * nquad0 + i];
1612 d0factors[3][j * nquad0 + i] +=
1613 df[3 * n][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1614 normal_3[0][j * nquad0 + i];
1615 d1factors[3][j * nquad0 + i] +=
1616 df[3 * n + 1][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1617 normal_3[0][j * nquad0 + i];
1618 d2factors[3][j * nquad0 + i] +=
1619 df[3 * n + 2][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1620 normal_3[0][j * nquad0 + i];
1626 for (
int j = 0; j < nquad2; ++j)
1628 for (
int i = 0; i < nquad1; ++i)
1630 d0factors[2][j * nquad1 + i] =
1631 df[0][j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1632 normal_2[0][j * nquad1 + i];
1633 d1factors[2][j * nquad1 + i] =
1634 df[1][j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1635 normal_2[0][j * nquad1 + i];
1636 d2factors[2][j * nquad1 + i] =
1637 df[2][j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1638 normal_2[0][j * nquad1 + i];
1640 d0factors[4][j * nquad1 + i] =
1641 df[0][j * nquad0 * nquad1 + i * nquad0] *
1642 normal_4[0][j * nquad1 + i];
1643 d1factors[4][j * nquad1 + i] =
1644 df[1][j * nquad0 * nquad1 + i * nquad0] *
1645 normal_4[0][j * nquad1 + i];
1646 d2factors[4][j * nquad1 + i] =
1647 df[2][j * nquad0 * nquad1 + i * nquad0] *
1648 normal_4[0][j * nquad1 + i];
1652 for (
int n = 1; n < ncoords; ++n)
1654 for (
int j = 0; j < nquad2; ++j)
1656 for (
int i = 0; i < nquad1; ++i)
1658 d0factors[2][j * nquad1 + i] +=
1659 df[3 * n][j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1660 normal_2[n][j * nquad1 + i];
1661 d1factors[2][j * nquad1 + i] +=
1663 [j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1664 normal_2[n][j * nquad1 + i];
1665 d2factors[2][j * nquad1 + i] +=
1667 [j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1668 normal_2[n][j * nquad1 + i];
1670 d0factors[4][j * nquad1 + i] +=
1671 df[3 * n][i * nquad0 + j * nquad0 * nquad1] *
1672 normal_4[n][j * nquad1 + i];
1673 d1factors[4][j * nquad1 + i] +=
1674 df[3 * n + 1][i * nquad0 + j * nquad0 * nquad1] *
1675 normal_4[n][j * nquad1 + i];
1676 d2factors[4][j * nquad1 + i] +=
1677 df[3 * n + 2][i * nquad0 + j * nquad0 * nquad1] *
1678 normal_4[n][j * nquad1 + i];
1686 for (
int i = 0; i < nquad0 * nquad1; ++i)
1688 d0factors[0][i] = df[0][0] * normal_0[0][i];
1689 d0factors[5][i] = df[0][0] * normal_5[0][i];
1691 d1factors[0][i] = df[1][0] * normal_0[0][i];
1692 d1factors[5][i] = df[1][0] * normal_5[0][i];
1694 d2factors[0][i] = df[2][0] * normal_0[0][i];
1695 d2factors[5][i] = df[2][0] * normal_5[0][i];
1698 for (
int n = 1; n < ncoords; ++n)
1700 for (
int i = 0; i < nquad0 * nquad1; ++i)
1702 d0factors[0][i] += df[3 * n][0] * normal_0[n][i];
1703 d0factors[5][i] += df[3 * n][0] * normal_5[n][i];
1705 d1factors[0][i] += df[3 * n + 1][0] * normal_0[n][i];
1706 d1factors[5][i] += df[3 * n + 1][0] * normal_5[n][i];
1708 d2factors[0][i] += df[3 * n + 2][0] * normal_0[n][i];
1709 d2factors[5][i] += df[3 * n + 2][0] * normal_5[n][i];
1714 for (
int i = 0; i < nquad0 * nquad2; ++i)
1716 d0factors[1][i] = df[0][0] * normal_1[0][i];
1717 d0factors[3][i] = df[0][0] * normal_3[0][i];
1719 d1factors[1][i] = df[1][0] * normal_1[0][i];
1720 d1factors[3][i] = df[1][0] * normal_3[0][i];
1722 d2factors[1][i] = df[2][0] * normal_1[0][i];
1723 d2factors[3][i] = df[2][0] * normal_3[0][i];
1726 for (
int n = 1; n < ncoords; ++n)
1728 for (
int i = 0; i < nquad0 * nquad2; ++i)
1730 d0factors[1][i] += df[3 * n][0] * normal_1[n][i];
1731 d0factors[3][i] += df[3 * n][0] * normal_3[n][i];
1733 d1factors[1][i] += df[3 * n + 1][0] * normal_1[n][i];
1734 d1factors[3][i] += df[3 * n + 1][0] * normal_3[n][i];
1736 d2factors[1][i] += df[3 * n + 2][0] * normal_1[n][i];
1737 d2factors[3][i] += df[3 * n + 2][0] * normal_3[n][i];
1742 for (
int i = 0; i < nquad1 * nquad2; ++i)
1744 d0factors[2][i] = df[0][0] * normal_2[0][i];
1745 d0factors[4][i] = df[0][0] * normal_4[0][i];
1747 d1factors[2][i] = df[1][0] * normal_2[0][i];
1748 d1factors[4][i] = df[1][0] * normal_4[0][i];
1750 d2factors[2][i] = df[2][0] * normal_2[0][i];
1751 d2factors[4][i] = df[2][0] * normal_4[0][i];
1754 for (
int n = 1; n < ncoords; ++n)
1756 for (
int i = 0; i < nquad1 * nquad2; ++i)
1758 d0factors[2][i] += df[3 * n][0] * normal_2[n][i];
1759 d0factors[4][i] += df[3 * n][0] * normal_4[n][i];
1761 d1factors[2][i] += df[3 * n + 1][0] * normal_2[n][i];
1762 d1factors[4][i] += df[3 * n + 1][0] * normal_4[n][i];
1764 d2factors[2][i] += df[3 * n + 2][0] * normal_2[n][i];
1765 d2factors[4][i] += df[3 * n + 2][0] * normal_4[n][i];
#define ASSERTL0(condition, msg)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Describes the specification for a Basis.
int GetNumPoints() const
Return points order at which basis is defined.
PointsKey GetPointsKey() const
Return distribution of points.
Defines a specification for a set of points.
virtual DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey) override
std::map< int, NormalVector > m_traceNormals
std::map< int, Array< OneD, NekDouble > > m_elmtBndNormDirElmtLen
the element length in each element boundary(Vertex, edge or face) normal direction calculated based o...
SpatialDomains::GeometrySharedPtr GetGeom() const
SpatialDomains::GeometrySharedPtr m_geom
void ComputeLaplacianMetric()
void ComputeGmatcdotMF(const Array< TwoD, const NekDouble > &df, const Array< OneD, const NekDouble > &direction, Array< OneD, Array< OneD, NekDouble > > &dfdir)
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
void ComputeQuadratureMetric()
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
const NormalVector & GetTraceNormal(const int id)
virtual void v_WeakDirectionalDerivMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray) override
Integrate the physical point list inarray over region.
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
virtual void v_MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
virtual void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
virtual DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey) override
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdRegions::StdMatrixKey &mkey) override
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
virtual void v_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
void v_DropLocStaticCondMatrix(const MatrixKey &mkey) override
virtual DNekScalMatSharedPtr v_GetLocMatrix(const MatrixKey &mkey) override
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
Calculate the inner product of inarray with respect to the given basis B = base0 * base1 * base2.
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...
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Calculates the inner product .
virtual void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdRegions::StdMatrixKey &mkey) override
void v_ComputeTraceNormal(const int face) override
HexExp(const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc, const SpatialDomains::HexGeomSharedPtr &geom)
Constructor using BasisKey class for quadrature points and order definition.
virtual void v_ComputeLaplacianMetric() override
virtual void v_ExtractDataToCoeffs(const NekDouble *data, const std::vector< unsigned int > &nummodes, const int mode_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType) override
virtual void v_LaplacianMatrixOp_MatFree_Kernel(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp) override
virtual 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_DropLocMatrix(const MatrixKey &mkey) override
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
virtual void v_GetCoord(const Array< OneD, const NekDouble > &Lcoords, Array< OneD, NekDouble > &coords) override
Retrieves the physical coordinates of a given set of reference coordinates.
virtual NekDouble v_PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals) override
This function evaluates the expansion at a single (arbitrary) point of the domain.
virtual void v_IProductWRTDirectionalDerivBase_SumFac(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
virtual void v_AlignVectorToCollapsedDir(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
virtual void v_PhysDirectionalDeriv(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &out) override
Physical derivative along a direction vector.
virtual NekDouble v_StdPhysEvaluate(const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals) override
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix(const MatrixKey &mkey) override
virtual void v_GetTracePhysMap(const int face, Array< OneD, int > &outarray) override
virtual LibUtilities::ShapeType v_DetShapeType() const override
Return the region shape using the enum-list of ShapeType.
virtual void v_ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
virtual void v_WeakDerivMatrixOp(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
virtual StdRegions::StdExpansionSharedPtr v_GetStdExp(void) const override
virtual StdRegions::StdExpansionSharedPtr v_GetLinStdExp(void) const override
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Calculate the inner product of inarray with respect to the elements basis.
virtual void v_NormalTraceDerivFactors(Array< OneD, Array< OneD, NekDouble > > &factors, Array< OneD, Array< OneD, NekDouble > > &d0factors, Array< OneD, Array< OneD, NekDouble > > &d1factors) override
: This method gets all of the factors which are required as part of the Gradient Jump Penalty stabili...
virtual void v_MassLevelCurvatureMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
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.
const LibUtilities::PointsKeyVector GetPointsKeys() const
const LibUtilities::BasisKey GetTraceBasisKey(const int i, int k=-1) const
This function returns the basis key belonging to the i-th trace.
void IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
this function calculates the inner product of a given function f with the different modes of the expa...
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.
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)
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
Class representing a hexehedral element in reference space.
MatrixType GetMatrixType() const
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.
void Interp3D(const BasisKey &fbasis0, const BasisKey &fbasis1, const BasisKey &fbasis2, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, const BasisKey &tbasis2, Array< OneD, NekDouble > &to)
this function interpolates a 3D function evaluated at the quadrature points of the 3D basis,...
void InterpCoeff3D(const BasisKey &fbasis0, const BasisKey &fbasis1, const BasisKey &fbasis2, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, const BasisKey &tbasis2, Array< OneD, NekDouble > &to)
void Interp2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
this function interpolates a 2D function evaluated at the quadrature points of the 2D basis,...
std::vector< PointsKey > PointsKeyVector
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
@ eOrtho_A
Principle Orthogonal Functions .
@ eGLL_Lagrange
Lagrange for SEM basis .
@ eModified_A
Principle Modified Functions .
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
GeomType
Indicates the type of element geometry.
@ eRegular
Geometry is straight-sided with constant geometric factors.
@ eMovingRegular
Currently unused.
@ eDeformed
Geometry is curved or has non-constant factors.
std::shared_ptr< HexGeom > HexGeomSharedPtr
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
@ eInvLaplacianWithUnityMean
std::shared_ptr< StdHexExp > StdHexExpSharedPtr
The above copyright notice and this permission notice shall be included.
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
static Array< OneD, NekDouble > NullNekDouble1DArray
std::shared_ptr< DNekMat > DNekMatSharedPtr
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
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 Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*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 Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/x.
void Vdiv(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 Zero(int n, T *x, const int incx)
Zero vector.
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
void Vvtvvtp(int n, const T *v, int incv, const T *w, int incw, const T *x, int incx, const T *y, int incy, T *z, int incz)
vvtvvtp (vector times vector plus vector times vector):
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
scalarT< T > sqrt(scalarT< T > in)