47 [[maybe_unused]]
int numcoeffs,
60 int nquad0 =
m_base[0]->GetNumPoints();
61 int nquad1 =
m_base[1]->GetNumPoints();
63 if (outarray_d0.size() > 0)
66 if (inarray.data() == outarray_d0.data())
69 Vmath::Vcopy(nquad0 * nquad1, inarray.get(), 1, wsp.get(), 1);
71 &(D0->GetPtr())[0], nquad0, &wsp[0], nquad0, 0.0,
72 &outarray_d0[0], nquad0);
77 &(D0->GetPtr())[0], nquad0, &inarray[0], nquad0, 0.0,
78 &outarray_d0[0], nquad0);
82 if (outarray_d1.size() > 0)
85 if (inarray.data() == outarray_d1.data())
88 Vmath::Vcopy(nquad0 * nquad1, inarray.get(), 1, wsp.get(), 1);
89 Blas::Dgemm(
'N',
'T', nquad0, nquad1, nquad1, 1.0, &wsp[0], nquad0,
90 &(D1->GetPtr())[0], nquad1, 0.0, &outarray_d1[0],
95 Blas::Dgemm(
'N',
'T', nquad0, nquad1, nquad1, 1.0, &inarray[0],
96 nquad0, &(D1->GetPtr())[0], nquad1, 0.0,
97 &outarray_d1[0], nquad0);
114 const int nq0 =
m_base[0]->GetNumPoints();
115 const int nq1 =
m_base[1]->GetNumPoints();
118 for (
int i = 0; i < nq1; ++i)
120 wsp[i] = StdExpansion::BaryEvaluate<0>(coll[0], &physvals[0] + i * nq0);
123 return StdExpansion::BaryEvaluate<1>(coll[1], &wsp[0]);
132 int nq0 =
m_base[0]->GetNumPoints();
133 int nq1 =
m_base[1]->GetNumPoints();
137 for (i = 0; i < nq1; ++i)
140 Blas::Ddot(nq0, &(I[0]->GetPtr())[0], 1, &physvals[i * nq0], 1);
144 val =
Blas::Ddot(nq1, I[1]->GetPtr(), 1, wsp1, 1);
152 [[maybe_unused]] std::array<NekDouble, 3> &firstOrderDerivs)
167 int nquad0 =
m_base[0]->GetNumPoints();
168 int nquad1 =
m_base[1]->GetNumPoints();
172 for (i = 0; i < nquad1; ++i)
174 Vmath::Vmul(nquad0, &inarray[0] + i * nquad0, 1, w0.get(), 1,
175 &tmp[0] + i * nquad0, 1);
178 for (i = 0; i < nquad0; ++i)
180 Vmath::Vmul(nquad1, &tmp[0] + i, nquad0, w1.get(), 1, &tmp[0] + i,
193 bool doCheckCollDir0,
bool doCheckCollDir1)
196 doCheckCollDir0, doCheckCollDir1);
204 bool doCheckCollDir0,
bool doCheckCollDir1)
207 doCheckCollDir0, doCheckCollDir1);
212 ASSERTL1((dir == 0) || (dir == 1),
"Invalid direction.");
214 int nquad0 =
m_base[0]->GetNumPoints();
215 int nquad1 =
m_base[1]->GetNumPoints();
216 int nqtot = nquad0 * nquad1;
217 int nmodes0 =
m_base[0]->GetNumModes();
226 for (
int i = 0; i < nqtot; i++)
230 m_base[1]->GetBdata(), tmp1, tmp3,
236 (*mat)(j, i) = tmp3[j];
241 for (
int i = 0; i < nqtot; i++)
245 m_base[1]->GetDbdata(), tmp1, tmp3,
251 (*mat)(j, i) = tmp3[j];
273 int nquad0 =
m_base[0]->GetNumPoints();
274 int nquad1 =
m_base[1]->GetNumPoints();
275 int nqtot = nquad0 * nquad1;
276 int nmodes0 =
m_base[0]->GetNumModes();
277 int nmodes1 =
m_base[1]->GetNumModes();
279 max(max(max(nqtot,
m_ncoeffs), nquad1 * nmodes0), nquad0 * nmodes1);
288 if (!(
m_base[0]->Collocation() &&
m_base[1]->Collocation()))
320 int nquad0 =
m_base[0]->GetNumPoints();
321 int nquad1 =
m_base[1]->GetNumPoints();
322 int nqtot = nquad0 * nquad1;
323 int nmodes0 =
m_base[0]->GetNumModes();
324 int nmodes1 =
m_base[1]->GetNumModes();
326 max(max(max(nqtot,
m_ncoeffs), nquad1 * nmodes0), nquad0 * nmodes1);
337 if (!(
m_base[0]->Collocation() &&
m_base[1]->Collocation()))
371 [[maybe_unused]]
const unsigned int traceid,
375 "This method must be defined at the individual shape level");
386 [[maybe_unused]]
int Q)
394 dir = (eid == 0) ? 0 : 1;
401 int numModes =
m_base[dir]->GetNumModes();
404 P = (
P == -1) ? numModes :
P;
407 if (maparray.size() !=
P)
413 for (i = 0; i <
P; ++i)
418 if (signarray.size() !=
P)
424 std::fill(signarray.get(), signarray.get() +
P, 1);
429 for (i = numModes; i <
P; ++i)
432 maparray[i] = maparray[0];
442 std::swap(maparray[0], maparray[1]);
444 for (i = 3; i < std::min(
P, numModes); i += 2)
452 ASSERTL1(
P == numModes,
"Different trace space edge dimension "
453 "and element edge dimension not currently "
454 "possible for GLL-Lagrange bases");
456 std::reverse(maparray.get(), maparray.get() +
P);
460 ASSERTL0(
false,
"Mapping not defined for this type of basis");
475 if (maparray.size() != map2.size())
480 for (
int i = 0; i < map2.size(); ++i)
482 maparray[i] = map1[map2[i]];
#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...
Describes the specification for a Basis.
void v_GenStdMatBwdDeriv(const int dir, DNekMatSharedPtr &mat) 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)=0
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
void v_GetTraceToElementMap(const int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient=eForwards, int P=-1, int Q=-1) 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)=0
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_GetTraceCoeffMap(const unsigned int traceid, Array< OneD, unsigned int > &maparray) override
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)
NekDouble v_PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals) override
This function evaluates the expansion at a single (arbitrary) point of the domain.
void v_GetElmtTraceToTraceMap(const unsigned int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient, int P, int Q) override
Determine the mapping to re-orientate the coefficients along the element trace (assumed to align with...
void LaplacianMatrixOp_MatFree_GenericImpl(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.
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void HelmholtzMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void LaplacianMatrixOp_MatFree_Kernel(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
NekDouble GetConstFactor(const ConstFactorType &factor) const
bool ConstFactorExists(const ConstFactorType &factor) const
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...
@ eModified_B
Principle Modified Functions .
@ P
Monomial polynomials .
@ eGauss_Lagrange
Lagrange Polynomials using the Gauss points.
@ eGLL_Lagrange
Lagrange for SEM basis .
@ eModified_A
Principle Modified Functions .
static const NekDouble kNekZeroTol
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 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.
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)