35 #include <boost/core/ignore_unused.hpp>
49 StdQuadExp::StdQuadExp()
58 :
StdExpansion(Ba.GetNumModes() * Bb.GetNumModes(), 2, Ba, Bb),
100 boost::ignore_unused(out_d2);
122 ASSERTL1(
false,
"input dir is out of range");
133 boost::ignore_unused(out_d2);
153 if (
m_base[0]->Collocation() &&
m_base[1]->Collocation())
156 inarray, 1, outarray, 1);
168 m_base[1]->GetNumModes());
171 outarray, wsp,
true,
true);
189 bool doCheckCollDir0,
bool doCheckCollDir1)
191 int nquad0 =
m_base[0]->GetNumPoints();
192 int nquad1 =
m_base[1]->GetNumPoints();
193 int nmodes0 =
m_base[0]->GetNumModes();
194 int nmodes1 =
m_base[1]->GetNumModes();
196 bool colldir0 = doCheckCollDir0 ? (
m_base[0]->Collocation()) :
false;
197 bool colldir1 = doCheckCollDir1 ? (
m_base[1]->Collocation()) :
false;
199 if (colldir0 && colldir1)
205 Blas::Dgemm(
'N',
'T', nquad0, nquad1, nmodes1, 1.0, &inarray[0], nquad0,
206 base1.get(), nquad1, 0.0, &outarray[0], nquad0);
210 Blas::Dgemm(
'N',
'N', nquad0, nmodes1, nmodes0, 1.0, base0.get(),
211 nquad0, &inarray[0], nmodes0, 0.0, &outarray[0], nquad0);
215 ASSERTL1(wsp.size() >= nquad0 * nmodes1,
216 "Workspace size is not sufficient");
220 Blas::Dgemm(
'N',
'N', nquad0, nmodes1, nmodes0, 1.0, base0.get(),
221 nquad0, &inarray[0], nmodes0, 0.0, &wsp[0], nquad0);
222 Blas::Dgemm(
'N',
'T', nquad0, nquad1, nmodes1, 1.0, &wsp[0], nquad0,
223 base1.get(), nquad1, 0.0, &outarray[0], nquad0);
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);
268 for (i = 0; i < 4; i++)
274 for (i = 0; i < npoints[0]; i++)
276 physEdge[0][i] = inarray[i];
277 physEdge[2][i] = inarray[npoints[0] * npoints[1] - 1 - i];
280 for (i = 0; i < npoints[1]; i++)
282 physEdge[1][i] = inarray[npoints[0] - 1 + i * npoints[0]];
284 inarray[(npoints[1] - 1) * npoints[0] - i * npoints[0]];
289 m_base[0]->GetBasisKey()),
291 m_base[1]->GetBasisKey())};
297 for (i = 0; i < 4; i++)
299 segexp[i % 2]->FwdTransBndConstrained(physEdge[i], coeffEdge[i]);
302 for (j = 0; j < nmodes[i % 2]; j++)
305 outarray[mapArray[j]] =
sign * coeffEdge[i][j];
325 int nInteriorDofs =
m_ncoeffs - nBoundaryDofs;
332 for (i = 0; i < nInteriorDofs; i++)
334 rhs[i] = tmp1[mapArray[i]];
337 Blas::Dgemv(
'N', nInteriorDofs, nInteriorDofs, 1.0,
338 &(matsys->GetPtr())[0], nInteriorDofs, rhs.get(), 1, 0.0,
341 for (i = 0; i < nInteriorDofs; i++)
343 outarray[mapArray[i]] = result[i];
380 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)
406 m_base[1]->GetBdata(), tmp, outarray, wsp,
413 m_base[1]->GetBdata(), inarray, outarray,
429 ASSERTL0((dir == 0) || (dir == 1),
"input dir is out of range");
431 int nquad0 =
m_base[0]->GetNumPoints();
432 int nquad1 =
m_base[1]->GetNumPoints();
433 int nqtot = nquad0 * nquad1;
434 int order0 =
m_base[0]->GetNumModes();
445 m_base[1]->GetDbdata(), tmp, outarray, wsp,
451 m_base[1]->GetBdata(), tmp, outarray, wsp,
471 bool doCheckCollDir0,
bool doCheckCollDir1)
473 int nquad0 =
m_base[0]->GetNumPoints();
474 int nquad1 =
m_base[1]->GetNumPoints();
475 int nmodes0 =
m_base[0]->GetNumModes();
476 int nmodes1 =
m_base[1]->GetNumModes();
478 bool colldir0 = doCheckCollDir0 ? (
m_base[0]->Collocation()) :
false;
479 bool colldir1 = doCheckCollDir1 ? (
m_base[1]->Collocation()) :
false;
481 if (colldir0 && colldir1)
487 Blas::Dgemm(
'N',
'N', nmodes0, nmodes1, nquad1, 1.0, inarray.get(),
488 nmodes0, base1.get(), nquad1, 0.0, outarray.get(), nmodes0);
492 Blas::Dgemm(
'T',
'N', nmodes0, nquad1, nquad0, 1.0, base0.get(), nquad0,
493 inarray.get(), nquad0, 0.0, outarray.get(), nmodes0);
497 ASSERTL1(wsp.size() >= nquad1 * nmodes0,
498 "Workspace size is not sufficient");
501 Blas::Dgemm(
'T',
'N', nmodes0, nquad1, nquad0, 1.0, base0.get(), nquad0,
502 inarray.get(), nquad0, 0.0, wsp.get(), nmodes0);
505 for (
int i = 0; i < nmodes0; ++i)
507 for (
int j = 0; j < nquad1; ++j)
509 wsp[j * nmodes0 + i] =
510 Blas::Ddot(nquad0, base0.get() + i * nquad0, 1,
511 inarray.get() + j * nquad0, 1);
515 Blas::Dgemm(
'N',
'N', nmodes0, nmodes1, nquad1, 1.0, wsp.get(), nmodes0,
516 base1.get(), nquad1, 0.0, outarray.get(), nmodes0);
547 int nquad0 =
m_base[0]->GetNumPoints();
548 int nquad1 =
m_base[1]->GetNumPoints();
551 int btmp0 =
m_base[0]->GetNumModes();
552 int mode0 = mode % btmp0;
553 int mode1 = mode / btmp0;
555 ASSERTL2(mode1 == (
int)floor((1.0 * mode) / btmp0),
556 "Integer Truncation not Equiv to Floor");
559 "calling argument mode is larger than total expansion order");
561 for (i = 0; i < nquad1; ++i)
564 &outarray[0] + i * nquad0, 1);
567 for (i = 0; i < nquad0; ++i)
570 &outarray[0] + i, nquad0, &outarray[0] + i, nquad0);
590 ASSERTL2((i >= 0) && (i <= 3),
"edge id is out of range");
592 if ((i == 0) || (i == 2))
604 ASSERTL2((i >= 0) && (i <= 4),
"edge id is out of range");
605 if ((i == 0) || (i == 2))
617 ASSERTL2((i >= 0) && (i <= 3),
"edge id is out of range");
619 if ((i == 0) || (i == 2))
632 boost::ignore_unused(j);
633 ASSERTL2((i >= 0) && (i <= 3),
"edge id is out of range");
635 if ((i == 0) || (i == 2))
655 "BasisType is not a boundary interior form");
659 "BasisType is not a boundary interior form");
669 "BasisType is not a boundary interior form");
673 "BasisType is not a boundary interior form");
679 const std::vector<unsigned int> &nummodes,
int &modes_offset)
681 int nmodes = nummodes[modes_offset] * nummodes[modes_offset + 1];
689 bool returnval =
false;
708 boost::ignore_unused(coords_2);
715 for (i = 0; i < nq1; ++i)
717 Blas::Dcopy(nq0, z0.get(), 1, &coords_0[0] + i * nq0, 1);
718 Vmath::Fill(nq0, z1[i], &coords_1[0] + i * nq0, 1);
742 const int nm0 =
m_base[0]->GetNumModes();
743 const int nm1 =
m_base[1]->GetNumModes();
745 return StdExpansion::BaryEvaluateBasis<0>(coords[0], mode % nm1) *
746 StdExpansion::BaryEvaluateBasis<1>(coords[1], mode / nm0);
757 int nummodes0, nummodes1;
758 int value1 = 0, value2 = 0;
764 nummodes0 =
m_base[0]->GetNumModes();
765 nummodes1 =
m_base[1]->GetNumModes();
777 value1 = 2 * nummodes0;
780 ASSERTL0(0,
"Mapping array is not defined for this expansion");
784 for (i = 0; i < value1; i++)
794 value2 = value1 + nummodes0 - 1;
800 ASSERTL0(0,
"Mapping array is not defined for this expansion");
804 for (i = 0; i < nummodes1 - 2; i++)
806 outarray[cnt++] = value1 + i * nummodes0;
807 outarray[cnt++] = value2 + i * nummodes0;
813 for (i = nummodes0 * (nummodes1 - 1); i <
GetNcoeffs(); i++)
824 int nummodes0, nummodes1;
831 nummodes0 =
m_base[0]->GetNumModes();
832 nummodes1 =
m_base[1]->GetNumModes();
840 startvalue = nummodes0;
843 startvalue = 2 * nummodes0;
846 ASSERTL0(0,
"Mapping array is not defined for this expansion");
859 ASSERTL0(0,
"Mapping array is not defined for this expansion");
863 for (i = 0; i < nummodes1 - 2; i++)
865 for (j = 0; j < nummodes0 - 2; j++)
867 outarray[cnt++] = startvalue + j;
869 startvalue += nummodes0;
877 if (useCoeffPacking ==
true)
879 switch (localVertexId)
890 localDOF =
m_base[0]->GetNumModes() - 1;
902 localDOF =
m_base[0]->GetNumModes() *
903 (
m_base[1]->GetNumModes() - 1);
907 localDOF =
m_base[0]->GetNumModes();
916 m_base[0]->GetNumModes() *
m_base[1]->GetNumModes() - 1;
920 localDOF =
m_base[0]->GetNumModes() + 1;
925 ASSERTL0(
false,
"eid must be between 0 and 3");
931 switch (localVertexId)
942 localDOF =
m_base[0]->GetNumModes() - 1;
955 m_base[0]->GetNumModes() *
m_base[1]->GetNumModes() - 1;
959 localDOF =
m_base[0]->GetNumModes() + 1;
967 localDOF =
m_base[0]->GetNumModes() *
968 (
m_base[1]->GetNumModes() - 1);
972 localDOF =
m_base[0]->GetNumModes();
977 ASSERTL0(
false,
"eid must be between 0 and 3");
991 ASSERTL1(traceid < 4,
"traceid must be between 0 and 3");
994 unsigned int order0 =
m_base[0]->GetNumModes();
995 unsigned int order1 =
m_base[1]->GetNumModes();
996 unsigned int numModes = (traceid % 2) ? order1 : order0;
998 if (maparray.size() != numModes)
1011 for (i = 0; i < numModes; i++)
1019 for (i = 0; i < numModes; i++)
1021 maparray[i] = i * order0 + 1;
1027 for (i = 0; i < numModes; i++)
1029 maparray[i] = order0 + i;
1035 for (i = 0; i < numModes; i++)
1037 maparray[i] = i * order0;
1052 for (i = 0; i < numModes; i++)
1060 for (i = 0; i < numModes; i++)
1062 maparray[i] = (i + 1) * order0 - 1;
1068 for (i = 0; i < numModes; i++)
1070 maparray[i] = order0 * (order1 - 1) + i;
1076 for (i = 0; i < numModes; i++)
1078 maparray[i] = order0 * i;
1088 ASSERTL0(
false,
"Mapping not defined for this type of basis");
1097 const int nummodes0 =
m_base[0]->GetNumModes();
1098 const int nummodes1 =
m_base[1]->GetNumModes();
1102 if (maparray.size() != nEdgeIntCoeffs)
1107 if (signarray.size() != nEdgeIntCoeffs)
1113 fill(signarray.get(), signarray.get() + nEdgeIntCoeffs, 1);
1122 for (i = 0; i < nEdgeIntCoeffs; i++)
1124 maparray[i] = i + 2;
1130 for (i = 0; i < nEdgeIntCoeffs; i++)
1132 maparray[i] = (i + 2) * nummodes0 + 1;
1138 for (i = 0; i < nEdgeIntCoeffs; i++)
1140 maparray[i] = nummodes0 + i + 2;
1146 for (i = 0; i < nEdgeIntCoeffs; i++)
1148 maparray[i] = (i + 2) * nummodes0;
1153 ASSERTL0(
false,
"eid must be between 0 and 3");
1159 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1171 for (i = 0; i < nEdgeIntCoeffs; i++)
1173 maparray[i] = i + 1;
1179 for (i = 0; i < nEdgeIntCoeffs; i++)
1181 maparray[i] = (i + 2) * nummodes0 - 1;
1187 for (i = 0; i < nEdgeIntCoeffs; i++)
1189 maparray[i] = nummodes0 * (nummodes1 - 1) + i + 1;
1195 for (i = 0; i < nEdgeIntCoeffs; i++)
1197 maparray[i] = nummodes0 * (i + 1);
1202 ASSERTL0(
false,
"eid must be between 0 and 3");
1207 reverse(maparray.get(), maparray.get() + nEdgeIntCoeffs);
1212 ASSERTL0(
false,
"Mapping not defined for this type of basis");
1233 int nq0 =
m_base[0]->GetNumPoints();
1234 int nq1 =
m_base[1]->GetNumPoints();
1257 for (
int i = 0; i < nq; ++i)
1259 for (
int j = 0; j < nq; ++j, ++cnt)
1262 coords[cnt][0] = -1.0 + 2 * j / (
NekDouble)(nq - 1);
1263 coords[cnt][1] = -1.0 + 2 * i / (
NekDouble)(nq - 1);
1267 for (
int i = 0; i < neq; ++i)
1271 I[0] =
m_base[0]->GetI(coll);
1272 I[1] =
m_base[1]->GetI(coll + 1);
1275 for (
int j = 0; j < nq1; ++j)
1281 Mat->GetRawPtr() + j * nq0 * neq + i, neq);
1293 for (i = 0; i < order1; ++i)
1295 (*Mat)(order0 * i + 1, i * order0 + 1) = 1.0;
1301 for (i = 0; i < order0; ++i)
1303 (*Mat)(order0 + i, order0 + i) = 1.0;
1317 (*Mat) = Imass * Iprod;
1325 int dir = (edge + 1) % 2;
1326 int nCoeffs =
m_base[dir]->GetNumModes();
1334 coords[0] = (edge == 0 || edge == 3) ? -1.0 : 1.0;
1341 Vmath::Vcopy(nCoeffs, m_Ix->GetPtr(), 1, Mat->GetPtr(), 1);
1367 int qa =
m_base[0]->GetNumPoints();
1368 int qb =
m_base[1]->GetNumPoints();
1369 int nmodes_a =
m_base[0]->GetNumModes();
1370 int nmodes_b =
m_base[1]->GetNumModes();
1371 int nmodes = min(nmodes_a, nmodes_b);
1384 OrthoExp.
FwdTrans(array, orthocoeffs);
1394 for (
int j = 0; j < nmodes_a; ++j)
1396 for (
int k = 0; k < nmodes_b; ++k)
1400 pow((1.0 * j) / (nmodes_a - 1), cutoff * nmodes_a),
1401 pow((1.0 * k) / (nmodes_b - 1), cutoff * nmodes_b));
1403 orthocoeffs[j * nmodes_b + k] *= SvvDiffCoeff * fac;
1414 max_ab = max(max_ab, 0);
1417 for (
int j = 0; j < nmodes_a; ++j)
1419 for (
int k = 0; k < nmodes_b; ++k)
1421 int maxjk = max(j, k);
1424 orthocoeffs[j * nmodes_b + k] *=
1434 min(nmodes_a, nmodes_b));
1440 for (
int j = 0; j < nmodes_a; ++j)
1442 for (
int k = 0; k < nmodes_b; ++k)
1444 if (j + k >= cutoff)
1446 orthocoeffs[j * nmodes_b + k] *=
1448 exp(-(j + k - nmodes) * (j + k - nmodes) /
1450 (j + k - cutoff + 1)))));
1454 orthocoeffs[j * nmodes_b + k] *= 0.0;
1462 OrthoExp.
BwdTrans(orthocoeffs, array);
1471 int qa =
m_base[0]->GetNumPoints();
1472 int qb =
m_base[1]->GetNumPoints();
1473 int nmodesA =
m_base[0]->GetNumModes();
1474 int nmodesB =
m_base[1]->GetNumModes();
1475 int P = nmodesA - 1;
1476 int Q = nmodesB - 1;
1487 int Pcut = cutoff *
P;
1488 int Qcut = cutoff * Q;
1492 OrthoExp.
FwdTrans(array, orthocoeffs);
1496 for (
int i = 0; i < nmodesA; ++i)
1498 for (
int j = 0; j < nmodesB; ++j)
1501 if (i > Pcut || j > Qcut)
1505 fac = max(fac1, fac2);
1506 fac = pow(fac, exponent);
1507 orthocoeffs[i * nmodesB + j] *= exp(-alpha * fac);
1513 OrthoExp.
BwdTrans(orthocoeffs, array);
1520 int n_coeffs = inarray.size();
1527 int nmodes0 =
m_base[0]->GetNumModes();
1528 int nmodes1 =
m_base[1]->GetNumModes();
1529 int numMax = nmodes0;
1549 for (
int i = 0; i < numMin + 1; ++i)
1551 Vmath::Vcopy(numMin, tmp = coeff + cnt, 1, tmp2 = coeff_tmp + cnt, 1);
1600 int nquad0 =
m_base[0]->GetNumPoints();
1601 int nquad1 =
m_base[1]->GetNumPoints();
1607 for (i = 0; i < nquad1; ++i)
1609 Vmath::Vmul(nquad0, inarray.get() + i * nquad0, 1, w0.get(), 1,
1610 outarray.get() + i * nquad0, 1);
1613 for (i = 0; i < nquad0; ++i)
1615 Vmath::Vmul(nquad1, outarray.get() + i, nquad0, w1.get(), 1,
1616 outarray.get() + i, nquad0);
1623 boost::ignore_unused(standard);
1625 int np1 =
m_base[0]->GetNumPoints();
1626 int np2 =
m_base[1]->GetNumPoints();
1627 int np = max(np1, np2);
1634 for (
int i = 0; i < np - 1; ++i)
1637 for (
int j = 0; j < np - 1; ++j)
1639 conn[cnt++] = row + j;
1640 conn[cnt++] = row + j + 1;
1641 conn[cnt++] = rowp1 + j;
1643 conn[cnt++] = rowp1 + j + 1;
1644 conn[cnt++] = rowp1 + j;
1645 conn[cnt++] = row + j + 1;
#define ASSERTL0(condition, msg)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
#define sign(a, b)
return the sign(b)*a
Describes the specification for a 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.
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
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)
virtual 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.
const LibUtilities::BasisSharedPtr & GetBasis(int dir) const
This function gets the shared point to basis in the dir direction.
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)
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
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.
int GetTraceNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th trace.
DNekMatSharedPtr GenMatrix(const StdMatrixKey &mkey)
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.
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void MassMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
const ConstFactorMap & GetConstFactors() const
MatrixType GetMatrixType() const
NekDouble GetConstFactor(const ConstFactorType &factor) const
bool ConstFactorExists(const ConstFactorType &factor) const
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_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
virtual void v_ExponentialFilter(Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff) override
virtual 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
virtual LibUtilities::ShapeType v_DetShapeType() const final override
virtual void v_FwdTransBndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
virtual void v_LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta) override
virtual int v_NumDGBndryCoeffs() const final override
virtual int v_GetVertexMap(int localVertexId, bool useCoeffPacking=false) override
virtual int v_GetNverts() const final override
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
virtual bool v_IsBoundaryInteriorExpansion() const override
virtual void v_ReduceOrderCoeffs(int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
virtual int v_GetTraceNcoeffs(const int i) const final override
virtual void v_HelmholtzMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
virtual void v_SVVLaplacianFilter(Array< OneD, NekDouble > &array, const StdMatrixKey &mkey) override
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_0, Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2) override
virtual NekDouble v_PhysEvaluateBasis(const Array< OneD, const NekDouble > &coords, int mode) override
This function evaluates the basis function mode mode at a point coords of the domain.
virtual int v_CalcNumberOfCoefficients(const std::vector< unsigned int > &nummodes, int &modes_offset) 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.
virtual void v_GetInteriorMap(Array< OneD, unsigned int > &outarray) override
NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray) override
Integrates the specified function over the domain.
virtual int v_GetTraceIntNcoeffs(const int i) const final override
virtual int v_NumBndryCoeffs() const final override
virtual 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_BwdTrans(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 StdMatrixKey &mkey) override
virtual void v_LocCollapsedToLocCoord(const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi) override
virtual const LibUtilities::BasisKey v_GetTraceBasisKey(const int i, const int j) const final override
virtual void v_GetTraceCoeffMap(const unsigned int traceid, Array< OneD, unsigned int > &maparray) override
Get the map of the coefficient location to teh local trace coefficients.
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey) override
virtual 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
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
virtual 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
virtual int v_GetTraceNumPoints(const int i) const final override
void v_GetTraceInteriorToElementMap(const int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation edgeOrient=eForwards) 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.
virtual void v_GetBoundaryMap(Array< OneD, unsigned int > &outarray) override
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override
virtual ~StdQuadExp() override
Destructor.
virtual void v_LaplacianMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
virtual void v_FillMode(const int mode, Array< OneD, NekDouble > &array) override
Fill outarray with mode mode of expansion.
virtual int v_GetNtraces() const final 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 = A x where A[m x n].
static void Dcopy(const int &n, const double *x, const int &incx, double *y, const int &incy)
BLAS level 1: Copy x to y.
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...
int getNumberOfCoefficients(int Na, int Nb)
BasisManagerT & BasisManager(void)
std::shared_ptr< Basis > BasisSharedPtr
void InterpCoeff2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
@ 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 .
@ eFourier
Fourier Expansion .
static const NekDouble kNekZeroTol
@ eFactorSVVDGKerDiffCoeff
@ eFactorSVVPowerKerDiffCoeff
const int kSVVDGFiltermodesmin
const int kSVVDGFiltermodesmax
const NekDouble kSVVDGFilter[9][11]
@ ePhysInterpToEquiSpaced
std::map< ConstFactorType, NekDouble > ConstFactorMap
std::shared_ptr< StdSegExp > StdSegExpSharedPtr
The above copyright notice and this permission notice shall be included.
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 Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
void Zero(int n, T *x, const int incx)
Zero vector.
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
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.