44 Ba.GetNumModes(), Bb.GetNumModes()),
47 Ba.GetNumModes(), Bb.GetNumModes()),
51 "order in 'a' direction is higher than order "
61 int nquad1 =
m_base[1]->GetNumPoints();
70 case LibUtilities::eGaussRadauMAlpha1Beta0:
79 for (i = 0; i < nquad1; ++i)
81 w1_tmp[i] = 0.5 * (1 - z1[i]) * w1[i];
112 int nquad0 =
m_base[0]->GetNumPoints();
113 int nquad1 =
m_base[1]->GetNumPoints();
123 if (out_d0.size() > 0)
127 for (i = 0; i < nquad1; ++i)
129 Blas::Dscal(nquad0, wsp[i], &out_d0[0] + i * nquad0, 1);
133 if (out_d1.size() > 0)
139 for (i = 0; i < nquad1; ++i)
141 Vmath::Vvtvp(nquad0, &wsp[0], 1, &out_d0[0] + i * nquad0, 1,
142 &out_d1[0] + i * nquad0, 1,
143 &out_d1[0] + i * nquad0, 1);
147 else if (out_d1.size() > 0)
152 for (i = 0; i < nquad1; ++i)
154 Blas::Dscal(nquad0, wsp[i], &diff0[0] + i * nquad0, 1);
160 for (i = 0; i < nquad1; ++i)
162 Vmath::Vvtvp(nquad0, &wsp[0], 1, &diff0[0] + i * nquad0, 1,
163 &out_d1[0] + i * nquad0, 1, &out_d1[0] + i * nquad0,
187 ASSERTL1(
false,
"input dir is out of range");
227 m_base[0]->GetNumModes());
238 [[maybe_unused]]
bool doCheckCollDir0,
239 [[maybe_unused]]
bool doCheckCollDir1)
243 int nquad0 =
m_base[0]->GetNumPoints();
244 int nquad1 =
m_base[1]->GetNumPoints();
245 int nmodes0 =
m_base[0]->GetNumModes();
246 int nmodes1 =
m_base[1]->GetNumModes();
248 ASSERTL1(wsp.size() >= nquad1 * nmodes0,
249 "Workspace size is not sufficient");
252 "Basis[1] is not of general tensor type");
254 for (i = mode = 0; i < nmodes0; ++i)
256 Blas::Dgemv(
'N', nquad1, nmodes1 - i, 1.0, base1.get() + mode * nquad1,
257 nquad1, &inarray[0] + mode, 1, 0.0, &wsp[0] + i * nquad1,
265 Blas::Daxpy(nquad1, inarray[1], base1.get() + nquad1, 1,
266 &wsp[0] + nquad1, 1);
269 Blas::Dgemm(
'N',
'T', nquad0, nquad1, nmodes0, 1.0, base0.get(), nquad0,
270 &wsp[0], nquad1, 0.0, &outarray[0], nquad0);
286 out = (*matsys) * in;
294 int npoints[2] = {
m_base[0]->GetNumPoints(),
m_base[1]->GetNumPoints()};
295 int nmodes[2] = {
m_base[0]->GetNumModes(),
m_base[1]->GetNumModes()};
297 fill(outarray.get(), outarray.get() +
m_ncoeffs, 0.0);
301 for (i = 0; i < 3; i++)
307 for (i = 0; i < npoints[0]; i++)
309 physEdge[0][i] = inarray[i];
312 for (i = 0; i < npoints[1]; i++)
314 physEdge[1][i] = inarray[npoints[0] - 1 + i * npoints[0]];
316 inarray[(npoints[1] - 1) * npoints[0] - i * npoints[0]];
321 m_base[0]->GetBasisKey()),
323 m_base[1]->GetBasisKey())};
329 for (i = 0; i < 3; i++)
331 segexp[i != 0]->FwdTransBndConstrained(physEdge[i], coeffEdge[i]);
334 for (j = 0; j < nmodes[i != 0]; j++)
337 outarray[mapArray[j]] =
sign * coeffEdge[i][j];
357 int nInteriorDofs =
m_ncoeffs - nBoundaryDofs;
364 for (i = 0; i < nInteriorDofs; i++)
366 rhs[i] = tmp1[mapArray[i]];
369 Blas::Dgemv(
'N', nInteriorDofs, nInteriorDofs, 1.0, &(matsys->GetPtr())[0],
370 nInteriorDofs, rhs.get(), 1, 0.0, result.get(), 1);
372 for (i = 0; i < nInteriorDofs; i++)
374 outarray[mapArray[i]] = result[i];
430 int nquad0 =
m_base[0]->GetNumPoints();
431 int nquad1 =
m_base[1]->GetNumPoints();
432 int order0 =
m_base[0]->GetNumModes();
434 if (multiplybyweights)
442 m_base[1]->GetBdata(), tmp, outarray, wsp);
448 m_base[1]->GetBdata(), inarray, outarray,
458 [[maybe_unused]]
bool doCheckCollDir0,
459 [[maybe_unused]]
bool doCheckCollDir1)
463 int nquad0 =
m_base[0]->GetNumPoints();
464 int nquad1 =
m_base[1]->GetNumPoints();
465 int nmodes0 =
m_base[0]->GetNumModes();
466 int nmodes1 =
m_base[1]->GetNumModes();
468 ASSERTL1(wsp.size() >= nquad1 * nmodes0,
469 "Workspace size is not sufficient");
471 Blas::Dgemm(
'T',
'N', nquad1, nmodes0, nquad0, 1.0, inarray.get(), nquad0,
472 base0.get(), nquad0, 0.0, wsp.get(), nquad1);
475 for (mode = i = 0; i < nmodes0; ++i)
477 Blas::Dgemv(
'T', nquad1, nmodes1 - i, 1.0, base1.get() + mode * nquad1,
478 nquad1, wsp.get() + i * nquad1, 1, 0.0,
479 outarray.get() + mode, 1);
487 Blas::Ddot(nquad1, base1.get() + nquad1, 1, wsp.get() + nquad1, 1);
503 int nquad0 =
m_base[0]->GetNumPoints();
504 int nquad1 =
m_base[1]->GetNumPoints();
505 int nqtot = nquad0 * nquad1;
506 int nmodes0 =
m_base[0]->GetNumModes();
507 int wspsize = max(max(nqtot,
m_ncoeffs), nquad1 * nmodes0);
515 for (i = 0; i < nquad1; ++i)
517 gfac0[i] = 2.0 / (1 - z1[i]);
520 for (i = 0; i < nquad1; ++i)
522 Vmath::Smul(nquad0, gfac0[i], &inarray[0] + i * nquad0, 1,
523 &tmp0[0] + i * nquad0, 1);
533 m_base[1]->GetBdata(), tmp0, outarray,
542 for (i = 0; i < nquad0; ++i)
544 gfac0[i] = 0.5 * (1 + z0[i]);
547 for (i = 0; i < nquad1; ++i)
549 Vmath::Vmul(nquad0, &gfac0[0], 1, &tmp0[0] + i * nquad0, 1,
550 &tmp0[0] + i * nquad0, 1);
554 m_base[1]->GetBdata(), tmp0, tmp3,
559 m_base[1]->GetDbdata(), tmp0, outarray,
567 ASSERTL1(
false,
"input dir is out of range");
592 eta[0] = 2. * (1. + xi[0]) / d1 - 1.0;
599 xi[0] = (1.0 + eta[0]) * (1.0 - eta[1]) * 0.5 - 1.0;
606 int nquad0 =
m_base[0]->GetNumPoints();
607 int nquad1 =
m_base[1]->GetNumPoints();
608 int order0 =
m_base[0]->GetNumModes();
609 int order1 =
m_base[1]->GetNumModes();
615 "total expansion order");
618 for (i = 0; i < order0; ++i, m += order1 - i)
634 for (i = 0; i < nquad1; ++i)
637 &outarray[0] + i * nquad0, 1);
641 for (i = 0; i < nquad0; ++i)
644 &outarray[0] + i, nquad0, &outarray[0] + i, nquad0);
656 const int nm1 =
m_base[1]->GetNumModes();
657 const double c = 1 + 2 * nm1;
658 const int mode0 = floor(0.5 * (c -
sqrt(c * c - 8 * mode)));
663 return StdExpansion::BaryEvaluateBasis<1>(coll[1], 1);
667 return StdExpansion::BaryEvaluateBasis<0>(coll[0], mode0) *
668 StdExpansion::BaryEvaluateBasis<1>(coll[1], mode);
674 std::array<NekDouble, 3> &firstOrderDerivs)
681 if ((1 - coll[1]) < 1e-5)
685 PhysDeriv(inarray, EphysDeriv0, EphysDeriv1);
688 I[0] =
GetBase()[0]->GetI(coll);
689 I[1] =
GetBase()[1]->GetI(coll + 1);
703 temp = firstOrderDerivs[0];
706 firstOrderDerivs[0] = firstOrderDerivs[0] * fac0;
709 NekDouble fac1 = fac0 * (coll[0] + 1) / 2;
712 firstOrderDerivs[1] += fac1 * temp;
735 "BasisType is not a boundary interior form");
737 "BasisType is not a boundary interior form");
745 "BasisType is not a boundary interior form");
747 "BasisType is not a boundary interior form");
754 ASSERTL2(i >= 0 && i <= 2,
"edge id is out of range");
768 ASSERTL2(i >= 0 && i <= 2,
"edge id is out of range");
782 ASSERTL2((i >= 0) && (i <= 2),
"edge id is out of range");
795 const std::vector<unsigned int> &nummodes,
int &modes_offset)
798 nummodes[modes_offset], nummodes[modes_offset + 1]);
814 for (i = 0; i < nq1; ++i)
816 for (j = 0; j < nq0; ++j)
818 coords_0[i * nq0 + j] = (1 + z0[j]) * (1 - z1[i]) / 2.0 - 1.0;
820 Vmath::Fill(nq0, z1[i], &coords_1[0] + i * nq0, 1);
831 const int i, [[maybe_unused]]
const int j)
const
833 ASSERTL2(i >= 0 && i <= 2,
"edge id is out of range");
847 return m_base[dir]->GetBasisKey();
853 "Unexpected points distribution " +
856 " in StdTriExp::v_GetTraceBasisKey");
865 case LibUtilities::eGaussRadauMAlpha1Beta0:
875 m_base[dir]->GetNumModes(),
898 m_base[dir]->GetNumModes(),
912 m_base[dir]->GetNumModes(),
919 "Unexpected points distribution " +
922 " in StdTriExp::v_GetTraceBasisKey");
931 case LibUtilities::eGaussRadauMAlpha1Beta0:
941 m_base[dir]->GetNumModes(),
948 "Unexpected points distribution " +
951 " in StdTriExp::v_GetTraceBasisKey");
959 "Information not available to set edge key");
973 "Mapping not defined for this type of basis");
976 if (useCoeffPacking ==
true)
978 switch (localVertexId)
992 localDOF =
m_base[1]->GetNumModes();
997 ASSERTL0(
false,
"eid must be between 0 and 2");
1004 switch (localVertexId)
1013 localDOF =
m_base[1]->GetNumModes();
1023 ASSERTL0(
false,
"eid must be between 0 and 2");
1036 "Expansion not of a proper type");
1040 int nummodes0, nummodes1;
1047 nummodes0 =
m_base[0]->GetNumModes();
1048 nummodes1 =
m_base[1]->GetNumModes();
1050 startvalue = 2 * nummodes1;
1052 for (i = 0; i < nummodes0 - 2; i++)
1054 for (j = 0; j < nummodes1 - 3 - i; j++)
1056 outarray[cnt++] = startvalue + j;
1058 startvalue += nummodes1 - 2 - i;
1066 "Expansion not of expected type");
1069 int nummodes0, nummodes1;
1077 nummodes0 =
m_base[0]->GetNumModes();
1078 nummodes1 =
m_base[1]->GetNumModes();
1080 value = 2 * nummodes1 - 1;
1081 for (i = 0; i < value; i++)
1087 for (i = 0; i < nummodes0 - 2; i++)
1089 outarray[cnt++] = value;
1090 value += nummodes1 - 2 - i;
1100 "Mapping not defined for this type of basis");
1102 ASSERTL1(eid < 3,
"eid must be between 0 and 2");
1105 int order0 =
m_base[0]->GetNumModes();
1106 int order1 =
m_base[1]->GetNumModes();
1107 int numModes = (eid == 0) ? order0 : order1;
1109 if (maparray.size() != numModes)
1119 for (i = 0; i < numModes; cnt += order1 - i, ++i)
1127 maparray[0] = order1;
1129 for (i = 2; i < numModes; i++)
1131 maparray[i] = order1 - 1 + i;
1137 for (i = 0; i < numModes; i++)
1144 ASSERTL0(
false,
"eid must be between 0 and 2");
1155 "Mapping not defined for this type of basis");
1157 const int nummodes1 =
m_base[1]->GetNumModes();
1160 if (maparray.size() != nEdgeIntCoeffs)
1165 if (signarray.size() != nEdgeIntCoeffs)
1171 fill(signarray.get(), signarray.get() + nEdgeIntCoeffs, 1);
1178 int cnt = 2 * nummodes1 - 1;
1179 for (i = 0; i < nEdgeIntCoeffs; cnt += nummodes1 - 2 - i, ++i)
1187 for (i = 0; i < nEdgeIntCoeffs; i++)
1189 maparray[i] = nummodes1 + 1 + i;
1195 for (i = 0; i < nEdgeIntCoeffs; i++)
1197 maparray[i] = 2 + i;
1203 ASSERTL0(
false,
"eid must be between 0 and 2");
1210 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1234 nq0 =
m_base[0]->GetNumPoints();
1235 nq1 =
m_base[1]->GetNumPoints();
1256 for (
int i = 0; i < nq; ++i)
1258 for (
int j = 0; j < nq - i; ++j, ++cnt)
1261 coords[cnt][0] = -1.0 + 2 * j / (
NekDouble)(nq - 1);
1262 coords[cnt][1] = -1.0 + 2 * i / (
NekDouble)(nq - 1);
1266 for (
int i = 0; i < neq; ++i)
1270 I[0] =
m_base[0]->GetI(coll);
1271 I[1] =
m_base[1]->GetI(coll + 1);
1274 for (
int j = 0; j < nq1; ++j)
1280 Mat->GetRawPtr() + j * nq0 * neq + i, neq);
1344 int qa =
m_base[0]->GetNumPoints();
1345 int qb =
m_base[1]->GetNumPoints();
1346 int nmodes_a =
m_base[0]->GetNumModes();
1347 int nmodes_b =
m_base[1]->GetNumModes();
1360 OrthoExp.
FwdTrans(array, orthocoeffs);
1371 for (
int j = 0; j < nmodes_a; ++j)
1373 for (
int k = 0; k < nmodes_b - j; ++k, ++cnt)
1376 pow((1.0 * j) / (nmodes_a - 1), cutoff * nmodes_a),
1377 pow((1.0 * k) / (nmodes_b - 1), cutoff * nmodes_b));
1379 orthocoeffs[cnt] *= (SvvDiffCoeff * fac);
1390 max_ab = max(max_ab, 0);
1394 for (
int j = 0; j < nmodes_a; ++j)
1396 for (
int k = 0; k < nmodes_b - j; ++k, ++cnt)
1398 int maxjk = max(j, k);
1401 orthocoeffs[cnt] *= SvvDiffCoeff *
kSVVDGFilter[max_ab][maxjk];
1410 min(nmodes_a, nmodes_b));
1413 int nmodes = min(nmodes_a, nmodes_b);
1418 for (
int j = 0; j < nmodes_a; ++j)
1420 for (
int k = 0; k < nmodes_b - j; ++k)
1422 if (j + k >= cutoff)
1426 exp(-(j + k - nmodes) * (j + k - nmodes) /
1427 ((
NekDouble)((j + k - cutoff + epsilon) *
1428 (j + k - cutoff + epsilon)))));
1432 orthocoeffs[cnt] *= 0.0;
1440 OrthoExp.
BwdTrans(orthocoeffs, array);
1447 int n_coeffs = inarray.size();
1448 int nquad0 =
m_base[0]->GetNumPoints();
1449 int nquad1 =
m_base[1]->GetNumPoints();
1454 int nqtot = nquad0 * nquad1;
1457 int nmodes0 =
m_base[0]->GetNumModes();
1458 int nmodes1 =
m_base[1]->GetNumModes();
1459 int numMin2 = nmodes0;
1480 m_TriExp->BwdTrans(inarray, phys_tmp);
1481 m_OrthoTriExp->FwdTrans(phys_tmp, coeff);
1483 for (i = 0; i < n_coeffs; i++)
1488 numMin += numMin2 - 1;
1493 m_OrthoTriExp->BwdTrans(coeff, phys_tmp);
1494 m_TriExp->FwdTrans(phys_tmp, outarray);
1506 int nquad0 =
m_base[0]->GetNumPoints();
1507 int nquad1 =
m_base[1]->GetNumPoints();
1514 for (i = 0; i < nquad1; ++i)
1516 Vmath::Vmul(nquad0, inarray.get() + i * nquad0, 1, w0.get(), 1,
1517 outarray.get() + i * nquad0, 1);
1523 case LibUtilities::eGaussRadauMAlpha1Beta0:
1524 for (i = 0; i < nquad1; ++i)
1526 Blas::Dscal(nquad0, 0.5 * w1[i], outarray.get() + i * nquad0,
1532 for (i = 0; i < nquad1; ++i)
1535 outarray.get() + i * nquad0, 1);
1544 int np1 =
m_base[0]->GetNumPoints();
1545 int np2 =
m_base[1]->GetNumPoints();
1546 int np = max(np1, np2);
1553 for (
int i = 0; i < np - 1; ++i)
1556 for (
int j = 0; j < np - i - 2; ++j)
1558 conn[cnt++] = row + j;
1559 conn[cnt++] = row + j + 1;
1560 conn[cnt++] = rowp1 + j;
1562 conn[cnt++] = rowp1 + j + 1;
1563 conn[cnt++] = rowp1 + j;
1564 conn[cnt++] = row + j + 1;
1567 conn[cnt++] = row + np - i - 2;
1568 conn[cnt++] = row + np - i - 1;
1569 conn[cnt++] = rowp1 + np - i - 2;
#define ASSERTL0(condition, msg)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
#define sign(a, b)
return the sign(b)*a
Describes the specification for a Basis.
int GetNumModes() const
Returns the order of the basis.
Defines a specification for a set of points.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d0, Array< OneD, NekDouble > &outarray_d1)
Calculate the 2D derivative in the local tensor/collapsed coordinate at the physical points.
void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
NekDouble BaryTensorDeriv(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
NekDouble Integral(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &w0, const Array< OneD, const NekDouble > &w1)
void BwdTrans_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)
void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
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)
The base class for all shapes.
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
void WeakDerivMatrixOp_MatFree(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
int NumBndryCoeffs(void) const
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
const Array< OneD, const LibUtilities::BasisSharedPtr > & GetBase() const
This function gets the shared point to basis.
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
void LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
NekDouble PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
This function evaluates the expansion at a single (arbitrary) point of the domain.
void 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 BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Backward transformation from coefficient space to physical space.
int GetTraceNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th trace.
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Forward transformation from physical space to coefficient space.
LibUtilities::NekManager< StdMatrixKey, DNekBlkMat, StdMatrixKey::opLess > m_stdStaticCondMatrixManager
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
void PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void MassMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
MatrixType GetMatrixType() const
NekDouble GetConstFactor(const ConstFactorType &factor) const
bool ConstFactorExists(const ConstFactorType &factor) const
void v_WeakDerivMatrixOp(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
LibUtilities::ShapeType v_DetShapeType() const final
int v_GetTraceNumPoints(const int i) const override
void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Backward tranform for triangular elements.
void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
NekDouble v_PhysEvaluate(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
void v_LocCollapsedToLocCoord(const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi) override
void v_GetCoords(Array< OneD, NekDouble > &coords_x, Array< OneD, NekDouble > &coords_y, Array< OneD, NekDouble > &coords_z) override
int v_CalcNumberOfCoefficients(const std::vector< unsigned int > &nummodes, int &modes_offset) override
void v_GetSimplexEquiSpacedConnectivity(Array< OneD, int > &conn, bool standard=true) override
void v_MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray) override
int v_GetNtraces() const final
void v_StdPhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray) override
void v_GetTraceInteriorToElementMap(const int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation edgeOrient=eForwards) override
int v_GetVertexMap(int localVertexId, bool useCoeffPacking=false) override
const LibUtilities::BasisKey v_GetTraceBasisKey(const int i, const int j) const override
DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey) override
void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) 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[p]*base1[pq] and put into ou...
int v_GetTraceNcoeffs(const int i) const override
void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray) override
Calculate the derivative of the physical points.
void v_GetTraceCoeffMap(const unsigned int traceid, Array< OneD, unsigned int > &maparray) override
void v_BwdTrans_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, bool doCheckCollDir1) override
bool v_IsBoundaryInteriorExpansion() const override
void v_FillMode(const int mode, Array< OneD, NekDouble > &outarray) override
void v_ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_FwdTransBndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
int v_GetTraceIntNcoeffs(const int i) const override
int v_NumDGBndryCoeffs() const override
void v_LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta) override
NekDouble v_PhysEvaluateBasis(const Array< OneD, const NekDouble > &coords, int mode) final
void v_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray) override
Integrates the specified function over the domain.
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.
void v_GetInteriorMap(Array< OneD, unsigned int > &outarray) override
int v_GetNverts() const final
void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey) override
void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
int v_NumBndryCoeffs() const override
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
void v_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, bool doCheckCollDir1) override
static void Dgemv(const char &trans, const int &m, const int &n, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = alpha A x plus beta y where A[m x n].
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
static double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = .
static void Dgemm(const char &transa, const char &transb, const int &m, const int &n, const int &k, const double &alpha, const double *a, const int &lda, const double *b, const int &ldb, const double &beta, double *c, const int &ldc)
BLAS level 3: Matrix-matrix multiply C = A x B where op(A)[m x k], op(B)[k x n], C[m x n] DGEMM perfo...
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
int getNumberOfCoefficients(int Na)
int getNumberOfCoefficients(int Na, int Nb)
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.
const std::string kPointsTypeStr[]
@ eGaussRadauMLegendre
1D Gauss-Radau-Legendre quadrature points, pinned at x=-1
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
@ ePolyEvenlySpaced
1D Evenly-spaced points using Lagrange polynomial
@ eModified_B
Principle Modified Functions .
@ eOrtho_A
Principle Orthogonal Functions .
@ eGLL_Lagrange
Lagrange for SEM basis .
@ eOrtho_B
Principle Orthogonal Functions .
@ eModified_A
Principle Modified Functions .
static const NekDouble kNekZeroTol
@ eFactorSVVDGKerDiffCoeff
@ eFactorSVVPowerKerDiffCoeff
const int kSVVDGFiltermodesmin
const int kSVVDGFiltermodesmax
const NekDouble kSVVDGFilter[9][11]
@ ePhysInterpToEquiSpaced
std::shared_ptr< StdTriExp > StdTriExpSharedPtr
std::shared_ptr< StdSegExp > StdSegExpSharedPtr
static Array< OneD, NekDouble > NullNekDouble1DArray
std::shared_ptr< DNekMat > DNekMatSharedPtr
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
void 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 Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
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)