50 : StdExpansion(Ba.GetNumModes() * Bb.GetNumModes(), 2, Ba, Bb),
51 StdExpansion2D(Ba.GetNumModes() * Bb.GetNumModes(), Ba, Bb),
55 std::string(
"QuadExpMatrix")),
56 m_staticCondMatrixManager(
std::bind(&
Expansion::CreateStaticCondMatrix,
57 this,
std::placeholders::_1),
58 std::string(
"QuadExpStaticCondMatrix"))
63 : StdExpansion(T), StdExpansion2D(T), StdQuadExp(T),
Expansion(T),
65 m_staticCondMatrixManager(T.m_staticCondMatrixManager)
71 int nquad0 =
m_base[0]->GetNumPoints();
72 int nquad1 =
m_base[1]->GetNumPoints();
80 Vmath::Vmul(nquad0 * nquad1, jac, 1, inarray, 1, tmp, 1);
84 Vmath::Smul(nquad0 * nquad1, jac[0], inarray, 1, tmp, 1);
88 ival = StdQuadExp::v_Integral(tmp);
97 int nquad0 =
m_base[0]->GetNumPoints();
98 int nquad1 =
m_base[1]->GetNumPoints();
99 int nqtot = nquad0 * nquad1;
105 StdQuadExp::v_PhysDeriv(inarray, diff0, diff1);
112 Vmath::Vvtvp(nqtot, df[1], 1, diff1, 1, out_d0, 1, out_d0, 1);
118 Vmath::Vvtvp(nqtot, df[3], 1, diff1, 1, out_d1, 1, out_d1, 1);
124 Vmath::Vvtvp(nqtot, df[5], 1, diff1, 1, out_d2, 1, out_d2, 1);
175 ASSERTL1(
false,
"input dir is out of range");
185 int nquad0 =
m_base[0]->GetNumPoints();
186 int nquad1 =
m_base[1]->GetNumPoints();
187 int nqtot = nquad0 * nquad1;
195 StdQuadExp::v_PhysDeriv(inarray, diff0, diff1);
202 for (
int i = 0; i < 2; ++i)
205 for (
int k = 0; k < (
m_geom->GetCoordim()); ++k)
207 Vmath::Vvtvp(nqtot, &df[2 * k + i][0], 1, &direction[k * nqtot],
208 1, &tangmat[i][0], 1, &tangmat[i][0], 1);
215 Vmath::Vmul(nqtot, &tangmat[0][0], 1, &diff0[0], 1, &out[0], 1);
216 Vmath::Vvtvp(nqtot, &tangmat[1][0], 1, &diff1[0], 1, &out[0], 1,
230 if ((
m_base[0]->Collocation()) && (
m_base[1]->Collocation()))
246 out = (*matsys) * in;
254 if ((
m_base[0]->Collocation()) && (
m_base[1]->Collocation()))
261 int npoints[2] = {
m_base[0]->GetNumPoints(),
m_base[1]->GetNumPoints()};
262 int nmodes[2] = {
m_base[0]->GetNumModes(),
m_base[1]->GetNumModes()};
264 fill(outarray.get(), outarray.get() +
m_ncoeffs, 0.0);
266 if (nmodes[0] == 1 && nmodes[1] == 1)
268 outarray[0] = inarray[0];
275 for (i = 0; i < 4; i++)
282 for (i = 0; i < npoints[0]; i++)
284 physEdge[0][i] = inarray[i];
285 physEdge[2][i] = inarray[npoints[0] * (npoints[1] - 1) + i];
288 for (i = 0; i < npoints[1]; i++)
290 physEdge[1][i] = inarray[npoints[0] - 1 + i * npoints[0]];
291 physEdge[3][i] = inarray[i * npoints[0]];
294 for (i = 0; i < 4; i++)
298 reverse((physEdge[i]).get(),
299 (physEdge[i]).get() + npoints[i % 2]);
304 for (i = 0; i < 4; i++)
314 for (i = 0; i < 4; i++)
316 segexp[i % 2]->FwdTransBndConstrained(physEdge[i], coeffEdge[i]);
319 for (j = 0; j < nmodes[i % 2]; j++)
322 outarray[mapArray[j]] =
sign * coeffEdge[i][j];
327 int nInteriorDofs =
m_ncoeffs - nBoundaryDofs;
329 if (nInteriorDofs > 0)
353 for (i = 0; i < nInteriorDofs; i++)
355 rhs[i] = tmp1[mapArray[i]];
358 Blas::Dgemv(
'N', nInteriorDofs, nInteriorDofs, matsys->Scale(),
359 &((matsys->GetOwnedMatrix())->GetPtr())[0],
360 nInteriorDofs, rhs.get(), 1, 0.0, result.get(), 1);
362 for (i = 0; i < nInteriorDofs; i++)
364 outarray[mapArray[i]] = result[i];
373 if (
m_base[0]->Collocation() &&
m_base[1]->Collocation())
394 int nquad0 =
m_base[0]->GetNumPoints();
395 int nquad1 =
m_base[1]->GetNumPoints();
396 int order0 =
m_base[0]->GetNumModes();
398 if (multiplybyweights)
404 StdQuadExp::IProductWRTBase_SumFacKernel(
m_base[0]->GetBdata(),
405 m_base[1]->GetBdata(), tmp,
406 outarray, wsp,
true,
true);
412 StdQuadExp::IProductWRTBase_SumFacKernel(
m_base[0]->GetBdata(),
413 m_base[1]->GetBdata(), inarray,
414 outarray, wsp,
true,
true);
422 ASSERTL1((dir == 0) || (dir == 1) || (dir == 2),
"Invalid direction.");
424 "Invalid direction.");
426 int nquad0 =
m_base[0]->GetNumPoints();
427 int nquad1 =
m_base[1]->GetNumPoints();
428 int nqtot = nquad0 * nquad1;
429 int nmodes0 =
m_base[0]->GetNumModes();
446 tmp1, tmp3, tmp4,
false,
true);
448 tmp2, outarray, tmp4,
true,
false);
456 ASSERTL1((dir == 0) || (dir == 1) || (dir == 2),
"Invalid direction.");
458 "Invalid direction.");
460 int nquad0 =
m_base[0]->GetNumPoints();
461 int nquad1 =
m_base[1]->GetNumPoints();
462 int nqtot = nquad0 * nquad1;
463 int nmodes0 =
m_base[0]->GetNumModes();
475 Vmath::Vmul(nqtot, &df[2 * dir][0], 1, inarray.get(), 1, tmp1.get(), 1);
476 Vmath::Vmul(nqtot, &df[2 * dir + 1][0], 1, inarray.get(), 1, tmp2.get(),
481 Vmath::Smul(nqtot, df[2 * dir][0], inarray.get(), 1, tmp1.get(), 1);
482 Vmath::Smul(nqtot, df[2 * dir + 1][0], inarray.get(), 1, tmp2.get(), 1);
491 int nq =
m_base[0]->GetNumPoints() *
m_base[1]->GetNumPoints();
500 Vmath::Vvtvvtp(nq, &normals[0][0], 1, &Fx[0], 1, &normals[1][0], 1,
501 &Fy[0], 1, &Fn[0], 1);
502 Vmath::Vvtvp(nq, &normals[2][0], 1, &Fz[0], 1, &Fn[0], 1, &Fn[0], 1);
506 Vmath::Svtsvtp(nq, normals[0][0], &Fx[0], 1, normals[1][0], &Fy[0], 1,
508 Vmath::Svtvp(nq, normals[2][0], &Fz[0], 1, &Fn[0], 1, &Fn[0], 1);
530 m_base[0]->GetPointsKey());
532 m_base[1]->GetPointsKey());
550 ASSERTL1(Lcoords[0] >= -1.0 && Lcoords[1] <= 1.0 && Lcoords[1] >= -1.0 &&
552 "Local coordinates are not in region [-1,1]");
555 for (i = 0; i <
m_geom->GetCoordim(); ++i)
557 coords[i] =
m_geom->GetCoord(i, Lcoords);
571 return StdExpansion2D::v_PhysEvaluate(Lcoord, physvals);
580 m_geom->GetLocCoords(coord, Lcoord);
582 return StdExpansion2D::v_PhysEvaluate(Lcoord, physvals);
587 std::array<NekDouble, 3> &firstOrderDerivs)
591 m_geom->GetLocCoords(coord, Lcoord);
592 return StdQuadExp::v_PhysEvaluate(Lcoord, inarray, firstOrderDerivs);
600 int nquad0 =
m_base[0]->GetNumPoints();
601 int nquad1 =
m_base[1]->GetNumPoints();
610 Vmath::Vcopy(nquad0, &(inarray[0]), 1, &(outarray[0]), 1);
613 Vmath::Vcopy(nquad1, &(inarray[0]) + (nquad0 - 1), nquad0,
617 Vmath::Vcopy(nquad0, &(inarray[0]) + nquad0 * (nquad1 - 1), 1,
621 Vmath::Vcopy(nquad1, &(inarray[0]), nquad0, &(outarray[0]), 1);
624 ASSERTL0(
false,
"edge value (< 3) is out of range");
634 if (
m_base[edge % 2]->GetPointsKey() !=
635 EdgeExp->GetBasis(0)->GetPointsKey())
642 EdgeExp->GetBasis(0)->GetPointsKey(), outarray);
652 Vmath::Reverse(EdgeExp->GetNumPoints(0), &outarray[0], 1, &outarray[0],
662 int nq0 =
m_base[0]->GetNumPoints();
663 int nq1 =
m_base[1]->GetNumPoints();
677 for (i = 0; i < nq0; i++)
680 Blas::Ddot(nq1, mat_gauss->GetOwnedMatrix()->GetPtr().get(),
681 1, &inarray[i], nq0);
687 for (i = 0; i < nq1; i++)
690 Blas::Ddot(nq0, mat_gauss->GetOwnedMatrix()->GetPtr().get(),
691 1, &inarray[i * nq0], 1);
697 for (i = 0; i < nq0; i++)
700 Blas::Ddot(nq1, mat_gauss->GetOwnedMatrix()->GetPtr().get(),
701 1, &inarray[i], nq0);
707 for (i = 0; i < nq1; i++)
710 Blas::Ddot(nq0, mat_gauss->GetOwnedMatrix()->GetPtr().get(),
711 1, &inarray[i * nq0], 1);
716 ASSERTL0(
false,
"edge value (< 3) is out of range");
723 int nquad0 =
m_base[0]->GetNumPoints();
724 int nquad1 =
m_base[1]->GetNumPoints();
731 for (
int i = 0; i < nquad0; ++i)
738 for (
int i = 0; i < nquad1; ++i)
740 outarray[i] = (nquad0 - 1) + i * nquad0;
745 for (
int i = 0; i < nquad0; ++i)
747 outarray[i] = i + nquad0 * (nquad1 - 1);
752 for (
int i = 0; i < nquad1; ++i)
754 outarray[i] = i * nquad0;
758 ASSERTL0(
false,
"edge value (< 3) is out of range");
767 int nquad0 =
m_base[0]->GetNumPoints();
768 int nquad1 =
m_base[1]->GetNumPoints();
794 for (i = 0; i < nquad0; ++i)
797 j[i] *
sqrt(g1[i] * g1[i] + g3[i] * g3[i]);
801 Vmath::Vcopy(nquad1, &(df[0][0]) + (nquad0 - 1), nquad0,
804 Vmath::Vcopy(nquad1, &(df[2][0]) + (nquad0 - 1), nquad0,
810 for (i = 0; i < nquad1; ++i)
813 j[i] *
sqrt(g0[i] * g0[i] + g2[i] * g2[i]);
818 Vmath::Vcopy(nquad0, &(df[1][0]) + (nquad0 * (nquad1 - 1)),
821 Vmath::Vcopy(nquad0, &(df[3][0]) + (nquad0 * (nquad1 - 1)),
824 Vmath::Vcopy(nquad0, &(jac[0]) + (nquad0 * (nquad1 - 1)), 1,
827 for (i = 0; i < nquad0; ++i)
830 j[i] *
sqrt(g1[i] * g1[i] + g3[i] * g3[i]);
839 for (i = 0; i < nquad1; ++i)
842 j[i] *
sqrt(g0[i] * g0[i] + g2[i] * g2[i]);
846 ASSERTL0(
false,
"edge value (< 3) is out of range");
852 int nqtot = nquad0 * nquad1;
873 for (i = 0; i < nquad0; ++i)
875 outarray[i] =
sqrt(g1_edge[i] * g1_edge[i] +
876 g3_edge[i] * g3_edge[i]);
888 for (i = 0; i < nquad1; ++i)
890 outarray[i] =
sqrt(g0_edge[i] * g0_edge[i] +
891 g2_edge[i] * g2_edge[i]);
904 for (i = 0; i < nquad0; ++i)
906 outarray[i] =
sqrt(g1_edge[i] * g1_edge[i] +
907 g3_edge[i] * g3_edge[i]);
921 for (i = 0; i < nquad1; ++i)
923 outarray[i] =
sqrt(g0_edge[i] * g0_edge[i] +
924 g2_edge[i] * g2_edge[i]);
931 ASSERTL0(
false,
"edge value (< 3) is out of range");
943 for (i = 0; i < nquad0; ++i)
945 outarray[i] = jac[0] *
sqrt(df[1][0] * df[1][0] +
946 df[3][0] * df[3][0]);
950 for (i = 0; i < nquad1; ++i)
952 outarray[i] = jac[0] *
sqrt(df[0][0] * df[0][0] +
953 df[2][0] * df[2][0]);
957 for (i = 0; i < nquad0; ++i)
959 outarray[i] = jac[0] *
sqrt(df[1][0] * df[1][0] +
960 df[3][0] * df[3][0]);
964 for (i = 0; i < nquad1; ++i)
966 outarray[i] = jac[0] *
sqrt(df[0][0] * df[0][0] +
967 df[2][0] * df[2][0]);
971 ASSERTL0(
false,
"edge value (< 3) is out of range");
985 for (i = 0; i < ptsKeys.size(); ++i)
996 geomFactors->GetDerivFactors(ptsKeys);
1007 for (i = 0; i < vCoordDim; ++i)
1026 for (i = 0; i < vCoordDim; ++i)
1028 Vmath::Fill(nqe, -df[2 * i + 1][0], normal[i], 1);
1032 for (i = 0; i < vCoordDim; ++i)
1038 for (i = 0; i < vCoordDim; ++i)
1044 for (i = 0; i < vCoordDim; ++i)
1050 ASSERTL0(
false,
"edge is out of range (edge < 4)");
1055 for (i = 0; i < vCoordDim; ++i)
1057 fac += normal[i][0] * normal[i][0];
1059 fac = 1.0 /
sqrt(fac);
1063 for (i = 0; i < vCoordDim; ++i)
1065 Vmath::Smul(nqe, fac, normal[i], 1, normal[i], 1);
1072 int nquad0 = ptsKeys[0].GetNumPoints();
1073 int nquad1 = ptsKeys[1].GetNumPoints();
1091 for (j = 0; j < nquad0; ++j)
1093 edgejac[j] = jac[j];
1094 for (i = 0; i < vCoordDim; ++i)
1096 normals[i * nquad0 + j] =
1097 -df[2 * i + 1][j] * edgejac[j];
1100 from_key = ptsKeys[0];
1103 for (j = 0; j < nquad1; ++j)
1105 edgejac[j] = jac[nquad0 * j + nquad0 - 1];
1106 for (i = 0; i < vCoordDim; ++i)
1108 normals[i * nquad1 + j] =
1109 df[2 * i][nquad0 * j + nquad0 - 1] * edgejac[j];
1112 from_key = ptsKeys[1];
1115 for (j = 0; j < nquad0; ++j)
1117 edgejac[j] = jac[nquad0 * (nquad1 - 1) + j];
1118 for (i = 0; i < vCoordDim; ++i)
1120 normals[i * nquad0 + j] =
1121 (df[2 * i + 1][nquad0 * (nquad1 - 1) + j]) *
1125 from_key = ptsKeys[0];
1128 for (j = 0; j < nquad1; ++j)
1130 edgejac[j] = jac[nquad0 * j];
1131 for (i = 0; i < vCoordDim; ++i)
1133 normals[i * nquad1 + j] =
1134 -df[2 * i][nquad0 * j] * edgejac[j];
1137 from_key = ptsKeys[1];
1140 ASSERTL0(
false,
"edge is out of range (edge < 3)");
1145 int nqtot = nquad0 * nquad1;
1152 for (j = 0; j < nquad0; ++j)
1154 for (i = 0; i < vCoordDim; ++i)
1156 Vmath::Vmul(nqtot, &(df[2 * i + 1][0]), 1, &jac[0],
1157 1, &(tmp_gmat[0]), 1);
1160 normals[i * nquad0 + j] = -tmp_gmat_edge[j];
1163 from_key = ptsKeys[0];
1166 for (j = 0; j < nquad1; ++j)
1168 for (i = 0; i < vCoordDim; ++i)
1170 Vmath::Vmul(nqtot, &(df[2 * i][0]), 1, &jac[0], 1,
1174 normals[i * nquad1 + j] = tmp_gmat_edge[j];
1177 from_key = ptsKeys[1];
1180 for (j = 0; j < nquad0; ++j)
1182 for (i = 0; i < vCoordDim; ++i)
1184 Vmath::Vmul(nqtot, &(df[2 * i + 1][0]), 1, &jac[0],
1185 1, &(tmp_gmat[0]), 1);
1188 normals[i * nquad0 + j] = tmp_gmat_edge[j];
1191 from_key = ptsKeys[0];
1194 for (j = 0; j < nquad1; ++j)
1196 for (i = 0; i < vCoordDim; ++i)
1198 Vmath::Vmul(nqtot, &(df[2 * i][0]), 1, &jac[0], 1,
1202 normals[i * nquad1 + j] = -tmp_gmat_edge[j];
1205 from_key = ptsKeys[1];
1208 ASSERTL0(
false,
"edge is out of range (edge < 3)");
1224 Vmath::Vmul(nqe, work, 1, normal[i], 1, normal[i], 1);
1231 Vmath::Vvtvp(nqe, normal[i], 1, normal[i], 1, work, 1, work, 1);
1241 Vmath::Vmul(nqe, normal[i], 1, work, 1, normal[i], 1);
1246 for (i = 0; i < vCoordDim; ++i)
1257 const NekDouble *data,
const std::vector<unsigned int> &nummodes,
1259 std::vector<LibUtilities::BasisType> &fromType)
1261 int data_order0 = nummodes[mode_offset];
1262 int fillorder0 = std::min(
m_base[0]->GetNumModes(), data_order0);
1264 int data_order1 = nummodes[mode_offset + 1];
1265 int order1 =
m_base[1]->GetNumModes();
1266 int fillorder1 = min(order1, data_order1);
1277 m_base[0]->GetPointsKey()),
1279 m_base[1]->GetPointsKey()));
1281 m_base[1]->GetBasisKey());
1303 "Extraction routine not set up for this basis");
1306 for (i = 0; i < fillorder0; ++i)
1308 Vmath::Vcopy(fillorder1, data + cnt, 1, coeffs + cnt1, 1);
1342 ASSERTL0(
false,
"basis is either not set up or not hierarchicial");
1348 return m_geom->GetEorient(edge);
1366 returnval = StdQuadExp::v_GenMatrix(mkey);
1378 return tmp->GetStdMatrix(mkey);
1405 StdExpansion::MassMatrixOp_MatFree(inarray, outarray, mkey);
1420 StdExpansion::LaplacianMatrixOp_MatFree(k1, k2, inarray, outarray, mkey);
1428 StdExpansion::WeakDerivMatrixOp_MatFree(i, inarray, outarray, mkey);
1435 StdExpansion::WeakDirectionalDerivMatrixOp_MatFree(inarray, outarray, mkey);
1442 StdExpansion::MassLevelCurvatureMatrixOp_MatFree(inarray, outarray, mkey);
1456 int n_coeffs = inarray.size();
1462 int nmodes0 =
m_base[0]->GetNumModes();
1463 int nmodes1 =
m_base[1]->GetNumModes();
1464 int numMax = nmodes0;
1482 for (
int i = 0; i < numMin + 1; ++i)
1484 Vmath::Vcopy(numMin, tmp = coeff + cnt, 1, tmp2 = coeff_tmp + cnt, 1);
1501 int nquad0 =
m_base[0]->GetNumPoints();
1502 int nquad1 =
m_base[1]->GetNumPoints();
1503 int nqtot = nquad0 * nquad1;
1504 int nmodes0 =
m_base[0]->GetNumModes();
1505 int nmodes1 =
m_base[1]->GetNumModes();
1507 max(max(max(nqtot,
m_ncoeffs), nquad1 * nmodes0), nquad0 * nmodes1);
1509 ASSERTL1(wsp.size() >= 3 * wspsize,
"Workspace is of insufficient size.");
1527 StdExpansion2D::PhysTensorDeriv(inarray, wsp1, wsp2);
1533 Vmath::Vvtvvtp(nqtot, &metric00[0], 1, &wsp1[0], 1, &metric01[0], 1,
1534 &wsp2[0], 1, &wsp0[0], 1);
1535 Vmath::Vvtvvtp(nqtot, &metric01[0], 1, &wsp1[0], 1, &metric11[0], 1,
1536 &wsp2[0], 1, &wsp2[0], 1);
1558 const unsigned int dim = 2;
1566 for (
unsigned int i = 0; i < dim; ++i)
1568 for (
unsigned int j = i; j < dim; ++j)
1607 StdQuadExp::v_SVVLaplacianFilter(array, mkey);
1629 if (d0factors.size() != 4)
1635 if (d0factors[0].size() != nquad0)
1643 if (d0factors[1].size() != nquad1)
1661 int ncoords = normal_0.size();
1668 for (
int i = 0; i < nquad0; ++i)
1670 d0factors[0][i] = df[0][i] * normal_0[0][i];
1671 d0factors[2][i] = df[0][nquad0 * (nquad1 - 1) + i] * normal_2[0][i];
1673 d1factors[0][i] = df[1][i] * normal_0[0][i];
1674 d1factors[2][i] = df[1][nquad0 * (nquad1 - 1) + i] * normal_2[0][i];
1677 for (
int n = 1; n < ncoords; ++n)
1681 for (
int i = 0; i < nquad0; ++i)
1683 d0factors[0][i] += df[2 * n][i] * normal_0[n][i];
1685 df[2 * n][nquad0 * (nquad1 - 1) + i] * normal_2[n][i];
1687 d1factors[0][i] += df[2 * n + 1][i] * normal_0[n][i];
1689 df[2 * n + 1][nquad0 * (nquad1 - 1) + i] * normal_2[n][i];
1694 for (
int i = 0; i < nquad1; ++i)
1696 d0factors[1][i] = df[0][(i + 1) * nquad0 - 1] * normal_1[0][i];
1697 d0factors[3][i] = df[0][i * nquad0] * normal_3[0][i];
1699 d1factors[1][i] = df[1][(i + 1) * nquad0 - 1] * normal_1[0][i];
1700 d1factors[3][i] = df[1][i * nquad0] * normal_3[0][i];
1703 for (
int n = 1; n < ncoords; ++n)
1705 for (
int i = 0; i < nquad1; ++i)
1708 df[2 * n][(i + 1) * nquad0 - 1] * normal_1[n][i];
1709 d0factors[3][i] += df[2 * n][i * nquad0] * normal_3[n][i];
1712 df[2 * n + 1][(i + 1) * nquad0 - 1] * normal_1[n][i];
1713 d1factors[3][i] += df[2 * n + 1][i * nquad0] * normal_3[n][i];
1720 for (
int i = 0; i < nquad0; ++i)
1722 d1factors[0][i] = df[1][0] * normal_0[0][i];
1723 d1factors[2][i] = df[1][0] * normal_2[0][i];
1727 for (
int i = 0; i < nquad1; ++i)
1729 d0factors[1][i] = df[0][0] * normal_1[0][i];
1730 d0factors[3][i] = df[0][0] * normal_3[0][i];
1733 for (
int n = 1; n < ncoords; ++n)
1737 for (
int i = 0; i < nquad0; ++i)
1739 d1factors[0][i] += df[2 * n + 1][0] * normal_0[n][i];
1740 d1factors[2][i] += df[2 * n + 1][0] * normal_2[n][i];
1745 for (
int i = 0; i < nquad1; ++i)
1747 d0factors[1][i] += df[2 * n][0] * normal_1[n][i];
1748 d0factors[3][i] += df[2 * n][0] * normal_3[n][i];
1754 for (
int i = 0; i < nquad0; ++i)
1756 d0factors[0][i] = df[0][0] * normal_0[0][i];
1757 d0factors[2][i] = df[0][0] * normal_2[0][i];
1761 for (
int i = 0; i < nquad1; ++i)
1763 d1factors[1][i] = df[1][0] * normal_1[0][i];
1764 d1factors[3][i] = df[1][0] * normal_3[0][i];
1767 for (
int n = 1; n < ncoords; ++n)
1771 for (
int i = 0; i < nquad0; ++i)
1773 d0factors[0][i] += df[2 * n][0] * normal_0[n][i];
1774 d0factors[2][i] += df[2 * n][0] * normal_2[n][i];
1779 for (
int i = 0; i < nquad1; ++i)
1781 d1factors[1][i] += df[2 * n + 1][0] * normal_1[n][i];
1782 d1factors[3][i] += df[2 * n + 1][0] * normal_3[n][i];
#define ASSERTL0(condition, msg)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
#define sign(a, b)
return the sign(b)*a
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.
size_t GetNumPoints() const
DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey) override
SpatialDomains::Geometry2DSharedPtr GetGeom2D() const
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
ExpansionSharedPtr GetLeftAdjacentElementExp() const
void ComputeLaplacianMetric()
void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
int GetLeftAdjacentElementTrace() const
void ComputeQuadratureMetric()
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
StdRegions::Orientation GetTraceOrient(int trace)
const NormalVector & GetTraceNormal(const int id)
void v_MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
StdRegions::Orientation v_GetTraceOrient(int edge) override
DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix(const MatrixKey &mkey) override
void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Transform a given function from physical quadrature space to coefficient space.
DNekScalMatSharedPtr v_GetLocMatrix(const MatrixKey &mkey) override
StdRegions::StdExpansionSharedPtr v_GetLinStdExp(void) const override
void v_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
void v_GetCoord(const Array< OneD, const NekDouble > &Lcoords, Array< OneD, NekDouble > &coords) override
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 and put into outarray.
void v_ComputeLaplacianMetric() override
void v_LaplacianMatrixOp_MatFree_Kernel(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp) override
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
void v_AlignVectorToCollapsedDir(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
void GetEdgeInterpVals(const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray) override
void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
NekDouble v_PhysEvaluate(const Array< OneD, const NekDouble > &coord, const Array< OneD, const NekDouble > &physvals) override
This function evaluates the expansion at a single (arbitrary) point of the domain.
void v_GetTracePhysVals(const int edge, const StdRegions::StdExpansionSharedPtr &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient) override
void v_WeakDirectionalDerivMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) 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=NullNekDouble1DArray) override
Calculate the derivative of the physical points.
NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray) override
Integrates the specified function over the domain.
void v_DropLocMatrix(const MatrixKey &mkey) override
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
void v_ExtractDataToCoeffs(const NekDouble *data, const std::vector< unsigned int > &nummodes, const int mode_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType) override
void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
void v_ComputeTraceNormal(const int edge) override
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.
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
StdRegions::StdExpansionSharedPtr v_GetStdExp(void) const override
DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey) override
void v_DropLocStaticCondMatrix(const MatrixKey &mkey) override
void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_MassLevelCurvatureMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
DNekMatSharedPtr v_CreateStdMatrix(const StdRegions::StdMatrixKey &mkey) override
void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdRegions::StdMatrixKey &mkey) override
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 (GJP) s...
NekDouble v_StdPhysEvaluate(const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals) override
void v_FwdTransBndConstrained(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_GetTracePhysMap(const int edge, Array< OneD, int > &outarray) override
void v_GetTraceQFactors(const int edge, Array< OneD, NekDouble > &outarray) override
void v_WeakDerivMatrixOp(const int i, 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 > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
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.
void HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
int NumBndryCoeffs(void) const
const LibUtilities::PointsKeyVector GetPointsKeys() const
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
const LibUtilities::BasisKey GetTraceBasisKey(const int i, int k=-1) const
This function returns the basis key belonging to the i-th trace.
void IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
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...
void LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void GetTraceToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
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 IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Array< OneD, LibUtilities::BasisSharedPtr > m_base
MatrixType GetMatrixType() const
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 double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = .
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 Interp1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
this function interpolates a 1D function evaluated at the quadrature points of the basis fbasis0 to ...
void InterpCoeff2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, 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
@ eGaussGaussLegendre
1D Gauss-Gauss-Legendre quadrature points
@ eGauss_Lagrange
Lagrange Polynomials using the Gauss points.
@ eOrtho_A
Principle Orthogonal Functions .
@ eGLL_Lagrange
Lagrange for SEM basis .
@ eModified_A
Principle Modified Functions .
std::shared_ptr< SegExp > SegExpSharedPtr
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< Geometry2D > Geometry2DSharedPtr
std::shared_ptr< StdQuadExp > StdQuadExpSharedPtr
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
@ eInvLaplacianWithUnityMean
std::map< ConstFactorType, NekDouble > ConstFactorMap
StdRegions::ConstFactorMap factors
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 Svtsvtp(int n, const T alpha, const T *x, int incx, const T beta, const T *y, int incy, T *z, int incz)
Svtsvtp (scalar times vector plus scalar times vector):
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 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 Reverse(int n, const T *x, const int incx, T *y, const int incy)
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)
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
scalarT< T > sqrt(scalarT< T > in)