34#include <boost/core/ignore_unused.hpp>
53 Ba.GetNumModes(), (Bb.GetNumModes())),
56 Ba.GetNumModes(), (Bb.GetNumModes())),
60 std::bind(&
Expansion2D::CreateMatrix, this, std::placeholders::_1),
61 std::string(
"TriExpMatrix")),
62 m_staticCondMatrixManager(std::bind(&
Expansion::CreateStaticCondMatrix,
63 this, std::placeholders::_1),
64 std::string(
"TriExpStaticCondMatrix"))
69 : StdExpansion(T), StdExpansion2D(T), StdTriExp(T),
Expansion(T),
71 m_staticCondMatrixManager(T.m_staticCondMatrixManager)
77 int nquad0 =
m_base[0]->GetNumPoints();
78 int nquad1 =
m_base[1]->GetNumPoints();
86 Vmath::Vmul(nquad0 * nquad1, jac, 1, inarray, 1, tmp, 1);
90 Vmath::Smul(nquad0 * nquad1, jac[0], inarray, 1, tmp, 1);
94 ival = StdTriExp::v_Integral(tmp);
103 int nquad0 =
m_base[0]->GetNumPoints();
104 int nquad1 =
m_base[1]->GetNumPoints();
105 int nqtot = nquad0 * nquad1;
112 StdTriExp::v_PhysDeriv(inarray, diff0, diff1);
119 Vmath::Vvtvp(nqtot, df[1], 1, diff1, 1, out_d0, 1, out_d0, 1);
125 Vmath::Vvtvp(nqtot, df[3], 1, diff1, 1, out_d1, 1, out_d1, 1);
131 Vmath::Vvtvp(nqtot, df[5], 1, diff1, 1, out_d2, 1, out_d2, 1);
182 ASSERTL1(
false,
"input dir is out of range");
197 int nquad0 =
m_base[0]->GetNumPoints();
198 int nquad1 =
m_base[1]->GetNumPoints();
199 int nqtot = nquad0 * nquad1;
208 StdTriExp::v_PhysDeriv(inarray, diff0, diff1);
216 for (
int i = 0; i < 2; ++i)
219 for (
int k = 0; k < (
m_geom->GetCoordim()); ++k)
221 Vmath::Vvtvp(nqtot, &df[2 * k + i][0], 1, &direction[k * nqtot],
222 1, &tangmat[i][0], 1, &tangmat[i][0], 1);
227 Vmath::Vmul(nqtot, &tangmat[0][0], 1, &diff0[0], 1, &out[0], 1);
228 Vmath::Vvtvp(nqtot, &tangmat[1][0], 1, &diff1[0], 1, &out[0], 1,
235 for (
int i = 0; i < 2; ++i)
238 for (
int k = 0; k < (
m_geom->GetCoordim()); ++k)
240 Vmath::Svtvp(nqtot, df[2 * k + i][0], &direction[k * nqtot], 1,
241 &tangmat[i][0], 1, &tangmat[i][0], 1);
246 Vmath::Vmul(nqtot, &tangmat[0][0], 1, &diff0[0], 1, &out[0], 1);
248 Vmath::Vvtvp(nqtot, &tangmat[1][0], 1, &diff1[0], 1, &out[0], 1,
266 out = (*matsys) * in;
274 int npoints[2] = {
m_base[0]->GetNumPoints(),
m_base[1]->GetNumPoints()};
275 int nmodes[2] = {
m_base[0]->GetNumModes(),
m_base[1]->GetNumModes()};
277 fill(outarray.get(), outarray.get() +
m_ncoeffs, 0.0);
279 if (nmodes[0] == 1 && nmodes[1] == 1)
281 outarray[0] = inarray[0];
287 for (i = 0; i < 3; i++)
295 for (i = 0; i < npoints[0]; i++)
297 physEdge[0][i] = inarray[i];
301 for (i = 0; i < npoints[1]; i++)
303 physEdge[1][i] = inarray[npoints[0] - 1 + i * npoints[0]];
304 physEdge[2][i] = inarray[i * npoints[0]];
313 for (i = 1; i < 3; i++)
321 for (i = 1; i < 3; i++)
327 m_base[0]->GetPointsKey(), physEdge[i]);
329 npoints[1] = npoints[0];
339 for (i = 0; i < 3; i++)
341 segexp[i]->FwdTransBndConstrained(physEdge[i], coeffEdge[i]);
346 for (j = 0; j < nmodes[i != 0]; j++)
349 outarray[mapArray[j]] =
sign * coeffEdge[i][j];
354 int nInteriorDofs =
m_ncoeffs - nBoundaryDofs;
356 if (nInteriorDofs > 0)
380 for (i = 0; i < nInteriorDofs; i++)
382 rhs[i] = tmp1[mapArray[i]];
385 Blas::Dgemv(
'N', nInteriorDofs, nInteriorDofs, matsys->Scale(),
386 &((matsys->GetOwnedMatrix())->GetPtr())[0], nInteriorDofs,
387 rhs.get(), 1, 0.0, result.get(), 1);
389 for (i = 0; i < nInteriorDofs; i++)
391 outarray[mapArray[i]] = result[i];
413 int nquad0 =
m_base[0]->GetNumPoints();
414 int nquad1 =
m_base[1]->GetNumPoints();
415 int order0 =
m_base[0]->GetNumModes();
417 if (multiplybyweights)
424 m_base[1]->GetBdata(), tmp, outarray, wsp);
431 m_base[1]->GetBdata(), inarray, outarray,
440 int nquad0 =
m_base[0]->GetNumPoints();
441 int nquad1 =
m_base[1]->GetNumPoints();
442 int nqtot = nquad0 * nquad1;
443 int nmodes0 =
m_base[0]->GetNumModes();
444 int wspsize = max(max(nqtot,
m_ncoeffs), nquad1 * nmodes0);
463 tmp2, outarray, tmp0);
471 ASSERTL1((dir == 0) || (dir == 1) || (dir == 2),
"Invalid direction.");
473 "Invalid direction.");
475 int nquad0 =
m_base[0]->GetNumPoints();
476 int nquad1 =
m_base[1]->GetNumPoints();
477 int nqtot = nquad0 * nquad1;
478 int nmodes0 =
m_base[0]->GetNumModes();
479 int wspsize = max(max(nqtot,
m_ncoeffs), nquad1 * nmodes0);
496 for (
int i = 0; i < nquad1; ++i)
498 gfac0[i] = 2.0 / (1 - z1[i]);
500 for (
int i = 0; i < nquad0; ++i)
502 gfac1[i] = 0.5 * (1 + z0[i]);
505 for (
int i = 0; i < nquad1; ++i)
507 Vmath::Smul(nquad0, gfac0[i], &inarray[0] + i * nquad0, 1,
508 &tmp0[0] + i * nquad0, 1);
511 for (
int i = 0; i < nquad1; ++i)
513 Vmath::Vmul(nquad0, &gfac1[0], 1, &tmp0[0] + i * nquad0, 1,
514 &tmp1[0] + i * nquad0, 1);
519 Vmath::Vmul(nqtot, &df[2 * dir][0], 1, &tmp0[0], 1, &tmp0[0], 1);
520 Vmath::Vmul(nqtot, &df[2 * dir + 1][0], 1, &tmp1[0], 1, &tmp1[0], 1);
521 Vmath::Vmul(nqtot, &df[2 * dir + 1][0], 1, &inarray[0], 1, &tmp2[0], 1);
525 Vmath::Smul(nqtot, df[2 * dir][0], tmp0, 1, tmp0, 1);
526 Vmath::Smul(nqtot, df[2 * dir + 1][0], tmp1, 1, tmp1, 1);
527 Vmath::Smul(nqtot, df[2 * dir + 1][0], inarray, 1, tmp2, 1);
551 int nquad0 =
m_base[0]->GetNumPoints();
552 int nquad1 =
m_base[1]->GetNumPoints();
553 int nqtot = nquad0 * nquad1;
554 int nmodes0 =
m_base[0]->GetNumModes();
555 int wspsize = max(max(nqtot,
m_ncoeffs), nquad1 * nmodes0);
571 for (i = 0; i < nquad1; ++i)
573 gfac0[i] = 2.0 / (1 - z1[i]);
575 for (i = 0; i < nquad0; ++i)
577 gfac1[i] = 0.5 * (1 + z0[i]);
579 for (i = 0; i < nquad1; ++i)
581 Vmath::Smul(nquad0, gfac0[i], &inarray[0] + i * nquad0, 1,
582 &tmp0[0] + i * nquad0, 1);
584 for (i = 0; i < nquad1; ++i)
586 Vmath::Vmul(nquad0, &gfac1[0], 1, &tmp0[0] + i * nquad0, 1,
587 &tmp1[0] + i * nquad0, 1);
594 Vmath::Vmul(nqtot, &dfdir[0][0], 1, &tmp0[0], 1, &tmp0[0], 1);
595 Vmath::Vmul(nqtot, &dfdir[1][0], 1, &tmp1[0], 1, &tmp1[0], 1);
596 Vmath::Vmul(nqtot, &dfdir[1][0], 1, &inarray[0], 1, &tmp2[0], 1);
598 Vmath::Vadd(nqtot, &tmp0[0], 1, &tmp1[0], 1, &tmp1[0], 1);
606 tmp2, outarray, tmp0);
615 int nq =
m_base[0]->GetNumPoints() *
m_base[1]->GetNumPoints();
624 Vmath::Vvtvvtp(nq, &normals[0][0], 1, &Fx[0], 1, &normals[1][0], 1,
625 &Fy[0], 1, &Fn[0], 1);
626 Vmath::Vvtvp(nq, &normals[2][0], 1, &Fz[0], 1, &Fn[0], 1, &Fn[0], 1);
630 Vmath::Svtsvtp(nq, normals[0][0], &Fx[0], 1, normals[1][0], &Fy[0], 1,
632 Vmath::Svtvp(nq, normals[2][0], &Fz[0], 1, &Fn[0], 1, &Fn[0], 1);
655 m_base[0]->GetPointsKey());
657 m_base[1]->GetPointsKey());
668 ASSERTL1(Lcoords[0] >= -1.0 && Lcoords[1] <= 1.0 && Lcoords[1] >= -1.0 &&
670 "Local coordinates are not in region [-1,1]");
674 for (i = 0; i <
m_geom->GetCoordim(); ++i)
676 coords[i] =
m_geom->GetCoord(i, Lcoords);
697 return StdExpansion2D::v_PhysEvaluate(Lcoord, physvals);
706 m_geom->GetLocCoords(coord, Lcoord);
708 return StdExpansion2D::v_PhysEvaluate(Lcoord, physvals);
713 std::array<NekDouble, 3> &firstOrderDerivs)
717 m_geom->GetLocCoords(coord, Lcoord);
718 return StdTriExp::v_PhysEvaluate(Lcoord, inarray, firstOrderDerivs);
726 int nquad0 =
m_base[0]->GetNumPoints();
727 int nquad1 =
m_base[1]->GetNumPoints();
734 Vmath::Vcopy(nquad0, &(inarray[0]), 1, &(outarray[0]), 1);
738 Vmath::Vcopy(nquad1, &(inarray[0]) + (nquad0 - 1), nquad0,
743 Vmath::Vcopy(nquad1, &(inarray[0]), nquad0, &(outarray[0]), 1);
747 ASSERTL0(
false,
"edge value (< 3) is out of range");
751 ASSERTL1(EdgeExp->GetBasis(0)->GetPointsType() ==
753 "Edge expansion should be GLL");
756 if (
m_base[edge ? 1 : 0]->GetPointsKey() !=
757 EdgeExp->GetBasis(0)->GetPointsKey())
764 EdgeExp->GetBasis(0)->GetPointsKey(), outarray);
775 Vmath::Reverse(EdgeExp->GetNumPoints(0), &outarray[0], 1, &outarray[0],
783 boost::ignore_unused(edge, outarray);
784 ASSERTL0(
false,
"Routine not implemented for triangular elements");
789 int nquad0 =
m_base[0]->GetNumPoints();
790 int nquad1 =
m_base[1]->GetNumPoints();
797 for (
int i = 0; i < nquad0; ++i)
804 for (
int i = 0; i < nquad1; ++i)
806 outarray[i] = (nquad0 - 1) + i * nquad0;
811 for (
int i = 0; i < nquad1; ++i)
813 outarray[i] = i * nquad0;
817 ASSERTL0(
false,
"edge value (< 3) is out of range");
829 for (i = 0; i < ptsKeys.size(); ++i)
841 geomFactors->GetDerivFactors(ptsKeys);
852 for (i = 0; i < dim; ++i)
879 Vmath::Fill(nqe, df[2 * i + 1][0] + df[2 * i][0], normal[i],
890 ASSERTL0(
false,
"Edge is out of range (edge < 3)");
897 fac += normal[i][0] * normal[i][0];
899 fac = 1.0 /
sqrt(fac);
912 int nquad0 = ptsKeys[0].GetNumPoints();
913 int nquad1 = ptsKeys[1].GetNumPoints();
926 for (j = 0; j < nquad0; ++j)
931 normals[i * nquad0 + j] =
932 -df[2 * i + 1][j] * edgejac[j];
935 from_key = ptsKeys[0];
938 for (j = 0; j < nquad1; ++j)
940 edgejac[j] = jac[nquad0 * j + nquad0 - 1];
943 normals[i * nquad1 + j] =
944 (df[2 * i][nquad0 * j + nquad0 - 1] +
945 df[2 * i + 1][nquad0 * j + nquad0 - 1]) *
949 from_key = ptsKeys[1];
952 for (j = 0; j < nquad1; ++j)
954 edgejac[j] = jac[nquad0 * j];
957 normals[i * nquad1 + j] =
958 -df[2 * i][nquad0 * j] * edgejac[j];
961 from_key = ptsKeys[1];
964 ASSERTL0(
false,
"edge is out of range (edge < 3)");
979 Vmath::Vmul(nqe, work, 1, normal[i], 1, normal[i], 1);
986 Vmath::Vvtvp(nqe, normal[i], 1, normal[i], 1, work, 1, work, 1);
996 Vmath::Vmul(nqe, normal[i], 1, work, 1, normal[i], 1);
1013 const NekDouble *data,
const std::vector<unsigned int> &nummodes,
1014 const int mode_offset,
NekDouble *coeffs,
1015 std::vector<LibUtilities::BasisType> &fromType)
1017 boost::ignore_unused(fromType);
1019 int data_order0 = nummodes[mode_offset];
1020 int fillorder0 = min(
m_base[0]->GetNumModes(), data_order0);
1021 int data_order1 = nummodes[mode_offset + 1];
1022 int order1 =
m_base[1]->GetNumModes();
1023 int fillorder1 = min(order1, data_order1);
1036 "Extraction routine not set up for this basis");
1039 for (i = 0; i < fillorder0; ++i)
1041 Vmath::Vcopy(fillorder1 - i, &data[cnt], 1, &coeffs[cnt1], 1);
1042 cnt += data_order1 - i;
1048 ASSERTL0(
false,
"basis is either not set up or not hierarchicial");
1072 returnval = StdTriExp::v_GenMatrix(mkey);
1086 return tmp->GetStdMatrix(mkey);
1113 StdExpansion::MassMatrixOp_MatFree(inarray, outarray, mkey);
1128 StdExpansion::LaplacianMatrixOp_MatFree(k1, k2, inarray, outarray, mkey);
1136 StdExpansion::WeakDerivMatrixOp_MatFree(i, inarray, outarray, mkey);
1143 StdExpansion::WeakDirectionalDerivMatrixOp_MatFree(inarray, outarray, mkey);
1150 StdExpansion::MassLevelCurvatureMatrixOp_MatFree(inarray, outarray, mkey);
1169 int nquad0 =
m_base[0]->GetNumPoints();
1170 int nquad1 =
m_base[1]->GetNumPoints();
1171 int nqtot = nquad0 * nquad1;
1172 int nmodes0 =
m_base[0]->GetNumModes();
1173 int nmodes1 =
m_base[1]->GetNumModes();
1175 max(max(max(nqtot,
m_ncoeffs), nquad1 * nmodes0), nquad0 * nmodes1);
1177 ASSERTL1(wsp.size() >= 3 * wspsize,
"Workspace is of insufficient size.");
1195 StdExpansion2D::PhysTensorDeriv(inarray, wsp1, wsp2);
1201 Vmath::Vvtvvtp(nqtot, &metric00[0], 1, &wsp1[0], 1, &metric01[0], 1,
1202 &wsp2[0], 1, &wsp0[0], 1);
1203 Vmath::Vvtvvtp(nqtot, &metric01[0], 1, &wsp1[0], 1, &metric11[0], 1,
1204 &wsp2[0], 1, &wsp2[0], 1);
1226 const unsigned int dim = 2;
1235 for (i = 0; i < dim; ++i)
1237 for (j = i; j < dim; ++j)
1245 const unsigned int nquad0 =
m_base[0]->GetNumPoints();
1246 const unsigned int nquad1 =
m_base[1]->GetNumPoints();
1250 for (i = 0; i < nquad1; i++)
1252 Blas::Dscal(nquad0, 2.0 / (1 - z1[i]), &dEta_dXi[0][0] + i * nquad0, 1);
1253 Blas::Dscal(nquad0, 2.0 / (1 - z1[i]), &dEta_dXi[1][0] + i * nquad0, 1);
1255 for (i = 0; i < nquad0; i++)
1257 Blas::Dscal(nquad1, 0.5 * (1 + z0[i]), &dEta_dXi[1][0] + i, nquad0);
1264 Vmath::Smul(nqtot, df[0][0], &dEta_dXi[0][0], 1, &tmp[0], 1);
1265 Vmath::Svtvp(nqtot, df[1][0], &dEta_dXi[1][0], 1, &tmp[0], 1, &tmp[0],
1273 Vmath::Smul(nqtot, df[2][0], &dEta_dXi[0][0], 1, &tmp[0], 1);
1274 Vmath::Svtvp(nqtot, df[3][0], &dEta_dXi[1][0], 1, &tmp[0], 1, &tmp[0],
1286 Vmath::Smul(nqtot, df[4][0], &dEta_dXi[0][0], 1, &tmp[0], 1);
1287 Vmath::Svtvp(nqtot, df[5][0], &dEta_dXi[1][0], 1, &tmp[0], 1,
1298 NekDouble g2 = df[1][0] * df[1][0] + df[3][0] * df[3][0];
1301 g2 += df[5][0] * df[5][0];
1308 Vmath::Vmul(nqtot, &df[0][0], 1, &dEta_dXi[0][0], 1, &tmp[0], 1);
1309 Vmath::Vvtvp(nqtot, &df[1][0], 1, &dEta_dXi[1][0], 1, &tmp[0], 1,
1319 Vmath::Vmul(nqtot, &df[2][0], 1, &dEta_dXi[0][0], 1, &tmp[0], 1);
1320 Vmath::Vvtvp(nqtot, &df[3][0], 1, &dEta_dXi[1][0], 1, &tmp[0], 1,
1335 Vmath::Vmul(nqtot, &df[4][0], 1, &dEta_dXi[0][0], 1, &tmp[0], 1);
1336 Vmath::Vvtvp(nqtot, &df[5][0], 1, &dEta_dXi[1][0], 1, &tmp[0], 1,
1351 for (
unsigned int i = 0; i < dim; ++i)
1353 for (
unsigned int j = i; j < dim; ++j)
1376 int n_coeffs = inarray.size();
1377 int nquad0 =
m_base[0]->GetNumPoints();
1378 int nquad1 =
m_base[1]->GetNumPoints();
1379 int nqtot = nquad0 * nquad1;
1380 int nmodes0 =
m_base[0]->GetNumModes();
1381 int nmodes1 =
m_base[1]->GetNumModes();
1382 int numMin2 = nmodes0, i;
1411 m_TriExp->BwdTrans(inarray, phys_tmp);
1412 m_OrthoTriExp->FwdTrans(phys_tmp, coeff);
1414 for (i = 0; i < n_coeffs; i++)
1419 numMin += numMin2 - 1;
1424 m_OrthoTriExp->BwdTrans(coeff, phys_tmp);
1425 m_TriExp->FwdTrans(phys_tmp, outarray);
1449 StdTriExp::v_SVVLaplacianFilter(array, mkey);
1465 boost::ignore_unused(d2factors);
1472 if (d0factors.size() != 3)
1478 if (d0factors[0].size() != nquad0)
1484 if (d0factors[1].size() != nquad1)
1500 int ncoords = normal_0.size();
1506 for (
int i = 0; i < nquad0; ++i)
1508 d1factors[0][i] = df[1][i] * normal_0[0][i];
1512 for (
int i = 0; i < nquad1; ++i)
1514 d0factors[1][i] = df[0][(i + 1) * nquad0 - 1] * normal_1[0][i];
1515 d0factors[2][i] = df[0][i * nquad0] * normal_2[0][i];
1518 for (
int n = 1; n < ncoords; ++n)
1522 for (
int i = 0; i < nquad0; ++i)
1524 d1factors[0][i] += df[2 * n + 1][i] * normal_0[n][i];
1529 for (
int i = 0; i < nquad1; ++i)
1532 df[2 * n][(i + 1) * nquad0 - 1] * normal_1[n][i];
1533 d0factors[2][i] += df[2 * n][i * nquad0] * normal_2[n][i];
1539 for (
int i = 0; i < nquad0; ++i)
1541 d0factors[0][i] = df[0][i] * normal_0[0][i];
1545 for (
int i = 0; i < nquad1; ++i)
1547 d1factors[1][i] = df[1][(i + 1) * nquad0 - 1] * normal_1[0][i];
1548 d1factors[2][i] = df[1][i * nquad0] * normal_2[0][i];
1551 for (
int n = 1; n < ncoords; ++n)
1555 for (
int i = 0; i < nquad0; ++i)
1557 d0factors[0][i] += df[2 * n][i] * normal_0[n][i];
1562 for (
int i = 0; i < nquad1; ++i)
1565 df[2 * n + 1][(i + 1) * nquad0 - 1] * normal_1[n][i];
1566 d1factors[2][i] += df[2 * n + 1][i * nquad0] * normal_2[n][i];
1573 for (
int i = 0; i < nquad0; ++i)
1575 d1factors[0][i] = df[1][0] * normal_0[0][i];
1579 for (
int i = 0; i < nquad1; ++i)
1581 d0factors[1][i] = df[0][0] * normal_1[0][i];
1582 d0factors[2][i] = df[0][0] * normal_2[0][i];
1585 for (
int n = 1; n < ncoords; ++n)
1589 for (
int i = 0; i < nquad0; ++i)
1591 d1factors[0][i] += df[2 * n + 1][0] * normal_0[n][i];
1596 for (
int i = 0; i < nquad1; ++i)
1598 d0factors[1][i] += df[2 * n][0] * normal_1[n][i];
1599 d0factors[2][i] += df[2 * n][0] * normal_2[n][i];
1605 for (
int i = 0; i < nquad0; ++i)
1607 d0factors[0][i] = df[0][0] * normal_0[0][i];
1611 for (
int i = 0; i < nquad1; ++i)
1613 d1factors[1][i] = df[1][0] * normal_1[0][i];
1614 d1factors[2][i] = df[1][0] * normal_2[0][i];
1617 for (
int n = 1; n < ncoords; ++n)
1621 for (
int i = 0; i < nquad0; ++i)
1623 d0factors[0][i] += df[2 * n][0] * normal_0[n][i];
1628 for (
int i = 0; i < nquad1; ++i)
1630 d1factors[1][i] += df[2 * n + 1][0] * normal_1[n][i];
1631 d1factors[2][i] += df[2 * n + 1][0] * normal_2[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
virtual 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 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
int GetLeftAdjacentElementTrace() const
void ComputeQuadratureMetric()
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
StdRegions::Orientation GetTraceOrient(int trace)
const NormalVector & GetTraceNormal(const int id)
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
virtual StdRegions::StdExpansionSharedPtr v_GetStdExp(void) const override
virtual void v_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
virtual void v_WeakDerivMatrixOp(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
virtual 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.
virtual void v_IProductWRTDerivBase(const int dir, 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
void v_DropLocMatrix(const MatrixKey &mkey) override
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) 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_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix(const MatrixKey &mkey) override
virtual void v_GetTracePhysMap(const int edge, Array< OneD, int > &outarray) override
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdRegions::StdMatrixKey &mkey) override
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
virtual void v_IProductWRTDirectionalDerivBase(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, 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 StdRegions::StdExpansionSharedPtr v_GetLinStdExp(void) const override
virtual void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
virtual void v_GetCoord(const Array< OneD, const NekDouble > &Lcoords, Array< OneD, NekDouble > &coords) override
virtual DNekScalMatSharedPtr v_GetLocMatrix(const MatrixKey &mkey) override
virtual StdRegions::Orientation v_GetTraceOrient(int edge) override
virtual NekDouble v_StdPhysEvaluate(const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals) override
virtual void v_ComputeTraceNormal(const int edge) override
virtual void v_GetTracePhysVals(const int edge, const StdRegions::StdExpansionSharedPtr &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient) override
virtual DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey) override
void v_DropLocStaticCondMatrix(const MatrixKey &mkey) 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 basis B=base0[p]*base1[pq] and put into ou...
virtual void v_GetTraceQFactors(const int edge, Array< OneD, NekDouble > &outarray) override
virtual void v_WeakDirectionalDerivMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
virtual void v_FwdTransBndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
virtual void v_ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) 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_MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
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_IProductWRTDirectionalDerivBase_SumFac(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Directinoal Derivative in the modal space in the dir direction of varcoeffs.
virtual void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdRegions::StdMatrixKey &mkey) 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=NullNekDouble1DArray) override
Calculate the derivative of the physical points.
virtual 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
virtual 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.
virtual void v_MassLevelCurvatureMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
TriExp(const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const SpatialDomains::Geometry2DSharedPtr &geom)
Constructor using BasisKey class for quadrature points and order definition.
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray) override
Integrates the specified function over the domain.
virtual void v_ComputeLaplacianMetric() 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 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)
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
void IProductWRTDirectionalDerivBase_SumFac(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
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)
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
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 void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
int getNumberOfCoefficients(int Na)
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 ...
std::vector< PointsKey > PointsKeyVector
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
@ eModified_B
Principle Modified Functions .
@ eOrtho_A
Principle Orthogonal Functions .
@ eOrtho_B
Principle Orthogonal Functions .
@ 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< StdExpansion > StdExpansionSharedPtr
@ eInvLaplacianWithUnityMean
std::shared_ptr< StdTriExp > StdTriExpSharedPtr
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 Svtsvtp(int n, const T alpha, const T *x, int incx, const T beta, const T *y, int incy, T *z, int incz)
svtvvtp (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)