49 if (lhsPointsKey < rhsPointsKey)
53 if (lhsPointsKey != rhsPointsKey)
101 std::shared_ptr<Basis> returnval(
new Basis(bkey));
102 returnval->Initialize();
110 "Cannot call Basis initialisation with zero or negative order");
112 "zero or negative numbers of points");
145 std::shared_ptr<NekMatrix<NekDouble>> ftB(
244 for (
p = 0;
p < numModes; ++
p, mode += numPoints)
249 scal =
sqrt(0.5 * (2.0 *
p + 1.0));
250 for (i = 0; i < numPoints; ++i)
257 Blas::Dgemm(
'n',
'n', numPoints, numModes, numPoints, 1.0, D,
258 numPoints,
m_bdata.data(), numPoints, 0.0,
279 for (
size_t p = 0;
p < numModes; ++
p)
281 for (
size_t q = 0;
q < numModes -
p; ++
q, mode += numPoints)
285 for (
size_t j = 0; j < numPoints; ++j)
288 sqrt(
p +
q + 1.0) * pow(0.5 * (1.0 -
z[j]),
p);
294 Blas::Dgemm(
'n',
'n', numPoints, numModes * (numModes + 1) / 2,
295 numPoints, 1.0, D, numPoints,
m_bdata.data(), numPoints,
313 size_t P = numModes - 1, Q = numModes - 1, R = numModes - 1;
316 for (
size_t p = 0;
p <=
P; ++
p)
318 for (
size_t q = 0;
q <= Q -
p; ++
q)
320 for (
size_t r = 0; r <= R -
p -
q; ++r, mode += numPoints)
323 2 *
p + 2 *
q + 2.0, 0.0);
324 for (
size_t k = 0; k < numPoints; ++k)
327 mode[k] *= pow(0.5 * (1.0 -
z[k]),
p +
q);
330 mode[k] *=
sqrt(r +
p +
q + 1.5);
338 numModes * (numModes + 1) * (numModes + 2) / 6,
339 numPoints, 1.0, D, numPoints,
m_bdata.data(), numPoints,
365 size_t P = numModes - 1, Q = numModes - 1, R = numModes - 1;
368 for (
size_t p = 0;
p <=
P; ++
p)
370 for (
size_t q = 0;
q <= Q; ++
q)
372 for (
size_t r = 0; r <= R - std::max(
p,
q);
373 ++r, mode += numPoints)
382 size_t pq = std::max(
p +
q,
size_t(0));
386 for (
size_t k = 0; k < numPoints; ++k)
389 mode[k] *= pow(0.5 * (1.0 -
z[k]), pq);
392 mode[k] *=
sqrt(r + pq + 1.5);
400 numModes * (numModes + 1) * (numModes + 2) / 6,
401 numPoints, 1.0, D, numPoints,
m_bdata.data(), numPoints,
413 for (i = 0; i < numPoints; ++i)
416 m_bdata[numPoints + i] = 0.5 * (1 +
z[i]);
419 mode =
m_bdata.data() + 2 * numPoints;
421 for (
p = 2;
p < numModes; ++
p, mode += numPoints)
426 for (i = 0; i < numPoints; ++i)
433 Blas::Dgemm(
'n',
'n', numPoints, numModes, numPoints, 1.0, D,
434 numPoints,
m_bdata.data(), numPoints, 0.0,
460 for (i = 0; i < numPoints; ++i)
462 m_bdata[0 * numPoints + i] = 0.5 * (1 -
z[i]);
463 m_bdata[1 * numPoints + i] = 0.5 * (1 +
z[i]);
466 mode =
m_bdata.data() + 2 * numPoints;
468 for (
q = 2;
q < numModes; ++
q, mode += numPoints)
473 for (i = 0; i < numPoints; ++i)
480 for (i = 0; i < numPoints; ++i)
482 mode[i] = 0.5 * (1 -
z[i]);
487 for (
q = 2;
q < numModes; ++
q, mode += numPoints)
492 for (i = 0; i < numPoints; ++i)
500 one_p_z =
m_bdata.data() + numPoints;
502 for (
p = 2;
p < numModes; ++
p)
504 for (i = 0; i < numPoints; ++i)
506 mode[i] =
m_bdata[i] * one_m_z_pow[i];
512 for (
q = 1;
q < numModes -
p; ++
q, mode += numPoints)
517 for (i = 0; i < numPoints; ++i)
519 mode[i] *= one_m_z_pow[i] * one_p_z[i];
524 Blas::Dgemm(
'n',
'n', numPoints, numModes * (numModes + 1) / 2,
525 numPoints, 1.0, D, numPoints,
m_bdata.data(), numPoints,
569 for (
p = 0;
p < numModes; ++
p)
571 N = numPoints * (numModes -
p) * (numModes -
p + 1) / 2;
574 B_offset += numPoints * (numModes -
p);
580 numModes * (numModes + 1) * (numModes + 2) / 6,
581 numPoints, 1.0, D, numPoints,
m_bdata.data(), numPoints,
619 N = numPoints * (numModes) * (numModes + 1) / 2;
623 B_offset += numPoints * (numModes);
625 N = numPoints * (numModes - 1);
631 N = numPoints * (numModes - 1) * (numModes) / 2;
636 B_offset += numPoints * (numModes - 1);
640 mode =
m_bdata.data() + offset;
642 for (
p = 2;
p < numModes; ++
p)
645 N = numPoints * (numModes -
p);
652 one_p_z =
m_bdata.data() + numPoints;
654 for (
q = 2;
q < numModes; ++
q)
657 for (i = 0; i < numPoints; ++i)
670 for (
size_t r = 1; r < numModes - std::max(
p,
q); ++r)
673 r - 1, 2 *
p + 2 *
q - 3, 1.0);
675 for (i = 0; i < numPoints; ++i)
677 mode[i] *= one_m_z_pow[i] * one_p_z[i];
686 numModes * (numModes + 1) * (2 * numModes + 1) / 6,
687 numPoints, 1.0, D, numPoints,
m_bdata.data(), numPoints,
695 std::shared_ptr<Points<NekDouble>>
m_points =
699 for (
p = 0;
p < numModes; ++
p, mode += numPoints)
701 for (
q = 0;
q < numPoints; ++
q)
709 Blas::Dgemm(
'n',
'n', numPoints, numModes, numPoints, 1.0, D,
710 numPoints,
m_bdata.data(), numPoints, 0.0,
719 std::shared_ptr<Points<NekDouble>>
m_points =
723 for (
p = 0;
p < numModes; ++
p, mode += numPoints)
725 for (
q = 0;
q < numPoints; ++
q)
733 Blas::Dgemm(
'n',
'n', numPoints, numModes, numPoints, 1.0, D,
734 numPoints,
m_bdata.data(), numPoints, 0.0,
743 "Fourier modes should be a factor of 2");
745 for (i = 0; i < numPoints; ++i)
753 for (
p = 1;
p < numModes / 2; ++
p)
755 for (i = 0; i < numPoints; ++i)
757 m_bdata[2 *
p * numPoints + i] = cos(
p * M_PI * (
z[i] + 1));
758 m_bdata[(2 *
p + 1) * numPoints + i] =
759 -sin(
p * M_PI * (
z[i] + 1));
762 -
p * M_PI * sin(
p * M_PI * (
z[i] + 1));
764 -
p * M_PI * cos(
p * M_PI * (
z[i] + 1));
773 for (i = 0; i < numPoints; ++i)
775 m_bdata[i] = cos(M_PI * (
z[i] + 1));
776 m_bdata[numPoints + i] = -sin(M_PI * (
z[i] + 1));
778 m_dbdata[i] = -M_PI * sin(M_PI * (
z[i] + 1));
779 m_dbdata[numPoints + i] = -M_PI * cos(M_PI * (
z[i] + 1));
782 for (
p = 1;
p < numModes / 2; ++
p)
784 for (i = 0; i < numPoints; ++i)
786 m_bdata[2 *
p * numPoints + i] = 0.;
787 m_bdata[(2 *
p + 1) * numPoints + i] = 0.;
790 m_dbdata[(2 *
p + 1) * numPoints + i] = 0.;
816 for (
p = 0, scal = 1;
p < numModes; ++
p, mode += numPoints)
821 for (i = 0; i < numPoints; ++i)
826 scal *= 4 * (
p + 1) * (
p + 1) / (2 *
p + 2) / (2 *
p + 1);
830 Blas::Dgemm(
'n',
'n', numPoints, numModes, numPoints, 1.0, D,
831 numPoints,
m_bdata.data(), numPoints, 0.0,
840 for (
size_t p = 0;
p < numModes; ++
p, mode += numPoints)
842 for (
size_t i = 0; i < numPoints; ++i)
844 mode[i] = pow(
z[i],
p);
849 Blas::Dgemm(
'n',
'n', numPoints, numModes, numPoints, 1.0, D,
850 numPoints,
m_bdata.data(), numPoints, 0.0,
857 "not implemented at this time.");
866 bool returnval =
false;
874 case eGaussKronrodLegendre:
#define ASSERTL0(condition, msg)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
int GetNumModes() const
Return order of basis from the basis specification.
Array< OneD, NekDouble > m_dbdata
Derivative Basis definition.
Basis()=delete
Private default constructor.
static std::shared_ptr< Basis > Create(const BasisKey &bkey)
Returns a new instance of a Basis with given BasisKey.
PointsSharedPtr m_points
Set of points.
int GetTotNumPoints() const
Return total number of points from the basis specification.
BasisKey m_basisKey
Basis specification.
BasisType GetBasisType() const
Return the type of expansion basis.
int GetNumPoints() const
Return the number of points from the basis specification.
std::shared_ptr< NekMatrix< NekDouble > > CalculateInterpMatrix(const BasisKey &tbasis0)
Calculate the interpolation Matrix for coefficient from one base (m_basisKey) to another (tbasis0)
NekManager< BasisKey, NekMatrix< NekDouble >, BasisKey::opLess > m_InterpManager
void GenBasis()
Generate appropriate basis and their derivatives.
Array< OneD, NekDouble > m_bdata
Basis definition.
static bool initBasisManager
Describes the specification for a Basis.
int GetNumPoints() const
Return points order at which basis is defined.
BasisType GetBasisType() const
Return type of expansion basis.
PointsKey GetPointsKey() const
Return distribution of points.
size_t m_nummodes
Expansion order.
int GetTotNumPoints() const
int GetTotNumModes() const
int GetNumModes() const
Returns the order of the basis.
PointsType GetPointsType() const
Return type of quadrature.
BasisType m_basistype
Expansion type.
bool Collocation() const
Determine if basis has collocation properties.
bool ExactIprodInt() const
Determine if basis has exact integration for inner product.
bool RegisterGlobalCreator(const CreateFuncType &createFunc)
Register the Global Create Function. The return value is just to facilitate calling statically.
Defines a specification for a set of points.
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
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...
BasisManagerT & BasisManager(void)
bool operator==(const BasisKey &x, const BasisKey &y)
std::shared_ptr< Basis > BasisSharedPtr
const char *const BasisTypeMap[]
PointsManagerT & PointsManager(void)
bool operator<(const BasisKey &lhs, const BasisKey &rhs)
bool operator>(const BasisKey &lhs, const BasisKey &rhs)
bool operator!=(const BasisKey &x, const BasisKey &y)
std::ostream & operator<<(std::ostream &os, const BasisKey &rhs)
@ eGaussRadauMLegendre
1D Gauss-Radau-Legendre quadrature points, pinned at x=-1
@ eGaussRadauKronrodMLegendre
1D Gauss-Radau-Kronrod-Legendre quadrature points, pinned at x=-1
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
@ eGaussLobattoKronrodLegendre
1D Lobatto Kronrod quadrature points
@ eGaussGaussLegendre
1D Gauss-Gauss-Legendre quadrature points
@ eGaussRadauPLegendre
1D Gauss-Radau-Legendre quadrature points, pinned at x=1
@ eModified_B
Principle Modified Functions .
@ P
Monomial polynomials .
@ eGauss_Lagrange
Lagrange Polynomials using the Gauss points.
@ eOrtho_A
Principle Orthogonal Functions .
@ eModified_C
Principle Modified Functions .
@ eGLL_Lagrange
Lagrange for SEM basis .
@ eFourierSingleMode
Fourier ModifiedExpansion with just the first mode .
@ eChebyshev
Chebyshev Polynomials.
@ eOrtho_C
Principle Orthogonal Functions .
@ eModifiedPyr_C
Principle Modified Functions.
@ eOrtho_B
Principle Orthogonal Functions .
@ eModified_A
Principle Modified Functions .
@ eFourierHalfModeIm
Fourier Modified expansions with just the imaginary part of the first mode .
@ eFourierHalfModeRe
Fourier Modified expansions with just the real part of the first mode .
@ eOrthoPyr_C
Principle Orthogonal Functions .
@ eFourier
Fourier Expansion .
std::vector< double > w(NPUPPER)
std::vector< double > z(NPUPPER)
std::vector< double > q(NPUPPER *NPUPPER)
double hgj(const int i, const double z, const double *zgj, const int np, const double alpha, const double beta)
Compute the value of the i th Lagrangian interpolant through the np Gauss-Jacobi points zgj at the ar...
double hglj(const int i, const double z, const double *zglj, const int np, const double alpha, const double beta)
Compute the value of the i th Lagrangian interpolant through the np Gauss-Lobatto-Jacobi points zgrj ...
void jacobfd(const int np, const double *z, double *poly_in, double *polyd, const int n, const double alpha, const double beta)
Routine to calculate Jacobi polynomials, , and their first derivative, .
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
scalarT< T > sqrt(scalarT< T > in)
bool operator()(const BasisKey &lhs, const BasisKey &rhs) const