41 #include <boost/math/special_functions/gamma.hpp>
45 namespace LibUtilities
55 if (lhsPointsKey < rhsPointsKey)
59 if (lhsPointsKey != rhsPointsKey)
107 std::shared_ptr<Basis> returnval(
new Basis(bkey));
108 returnval->Initialize();
117 "Cannot call Basis initialisation with zero or negative order");
119 "zero or negative numbers of points");
152 std::shared_ptr<NekMatrix<NekDouble>> ftB(
235 D = &(
m_points->GetD()->GetPtr())[0];
251 for (
p = 0;
p < numModes; ++
p, mode += numPoints)
255 scal =
sqrt(0.5 * (2.0 *
p + 1.0));
256 for (i = 0; i < numPoints; ++i)
262 Blas::Dgemm(
'n',
'n', numPoints, numModes, numPoints, 1.0, D,
263 numPoints,
m_bdata.data(), numPoints, 0.0,
282 for (
int p = 0;
p < numModes; ++
p)
284 for (
int q = 0; q < numModes -
p; ++q, mode += numPoints)
288 for (
int j = 0; j < numPoints; ++j)
291 sqrt(
p + q + 1.0) * pow(0.5 * (1.0 - z[j]),
p);
297 Blas::Dgemm(
'n',
'n', numPoints, numModes * (numModes + 1) / 2,
298 numPoints, 1.0, D, numPoints,
m_bdata.data(), numPoints,
316 int P = numModes - 1, Q = numModes - 1, R = numModes - 1;
319 for (
int p = 0;
p <=
P; ++
p)
321 for (
int q = 0; q <= Q -
p; ++q)
323 for (
int r = 0; r <= R -
p - q; ++r, mode += numPoints)
326 2 *
p + 2 * q + 2.0, 0.0);
327 for (
int k = 0; k < numPoints; ++k)
330 mode[k] *= pow(0.5 * (1.0 - z[k]),
p + q);
333 mode[k] *=
sqrt(r +
p + q + 1.5);
341 numModes * (numModes + 1) * (numModes + 2) / 6,
342 numPoints, 1.0, D, numPoints,
m_bdata.data(), numPoints,
368 int P = numModes - 1, Q = numModes - 1, R = numModes - 1;
371 for (
int p = 0;
p <=
P; ++
p)
373 for (
int q = 0; q <= Q; ++q)
375 for (
int r = 0; r <= R - std::max(
p, q);
376 ++r, mode += numPoints)
385 int pq = std::max(
p + q, 0);
389 for (
int k = 0; k < numPoints; ++k)
392 mode[k] *= pow(0.5 * (1.0 - z[k]), pq);
395 mode[k] *=
sqrt(r + pq + 1.5);
403 numModes * (numModes + 1) * (numModes + 2) / 6,
404 numPoints, 1.0, D, numPoints,
m_bdata.data(), numPoints,
416 for (i = 0; i < numPoints; ++i)
419 m_bdata[numPoints + i] = 0.5 * (1 + z[i]);
422 mode =
m_bdata.data() + 2 * numPoints;
424 for (
p = 2;
p < numModes; ++
p, mode += numPoints)
429 for (i = 0; i < numPoints; ++i)
436 Blas::Dgemm(
'n',
'n', numPoints, numModes, numPoints, 1.0, D,
437 numPoints,
m_bdata.data(), numPoints, 0.0,
463 for (i = 0; i < numPoints; ++i)
465 m_bdata[0 * numPoints + i] = 0.5 * (1 - z[i]);
466 m_bdata[1 * numPoints + i] = 0.5 * (1 + z[i]);
469 mode =
m_bdata.data() + 2 * numPoints;
471 for (q = 2; q < numModes; ++q, mode += numPoints)
476 for (i = 0; i < numPoints; ++i)
483 for (i = 0; i < numPoints; ++i)
485 mode[i] = 0.5 * (1 - z[i]);
490 for (q = 2; q < numModes; ++q, mode += numPoints)
495 for (i = 0; i < numPoints; ++i)
503 one_p_z =
m_bdata.data() + numPoints;
505 for (
p = 2;
p < numModes; ++
p)
507 for (i = 0; i < numPoints; ++i)
509 mode[i] =
m_bdata[i] * one_m_z_pow[i];
515 for (q = 1; q < numModes -
p; ++q, mode += numPoints)
520 for (i = 0; i < numPoints; ++i)
522 mode[i] *= one_m_z_pow[i] * one_p_z[i];
527 Blas::Dgemm(
'n',
'n', numPoints, numModes * (numModes + 1) / 2,
528 numPoints, 1.0, D, numPoints,
m_bdata.data(), numPoints,
572 for (
p = 0;
p < numModes; ++
p)
574 N = numPoints * (numModes -
p) * (numModes -
p + 1) / 2;
577 B_offset += numPoints * (numModes -
p);
583 numModes * (numModes + 1) * (numModes + 2) / 6,
584 numPoints, 1.0, D, numPoints,
m_bdata.data(), numPoints,
622 N = numPoints * (numModes) * (numModes + 1) / 2;
626 B_offset += numPoints * (numModes);
628 N = numPoints * (numModes - 1);
634 N = numPoints * (numModes - 1) * (numModes) / 2;
639 B_offset += numPoints * (numModes - 1);
644 mode =
m_bdata.data() + offset;
646 for (
p = 2;
p < numModes; ++
p)
649 N = numPoints * (numModes -
p);
656 one_p_z =
m_bdata.data() + numPoints;
658 for (q = 2; q < numModes; ++q)
661 for (i = 0; i < numPoints; ++i)
667 mode[i] = pow(
m_bdata[i],
p + q - 2);
674 for (
int r = 1; r < numModes - std::max(
p, q); ++r)
677 2 *
p + 2 * q - 3, 1.0);
679 for (i = 0; i < numPoints; ++i)
681 mode[i] *= one_m_z_pow[i] * one_p_z[i];
690 numModes * (numModes + 1) * (2 * numModes + 1) / 6,
691 numPoints, 1.0, D, numPoints,
m_bdata.data(), numPoints,
699 std::shared_ptr<Points<NekDouble>>
m_points =
703 for (
p = 0;
p < numModes; ++
p, mode += numPoints)
705 for (q = 0; q < numPoints; ++q)
713 Blas::Dgemm(
'n',
'n', numPoints, numModes, numPoints, 1.0, D,
714 numPoints,
m_bdata.data(), numPoints, 0.0,
722 std::shared_ptr<Points<NekDouble>>
m_points =
726 for (
p = 0;
p < numModes; ++
p, mode += numPoints)
728 for (q = 0; q < numPoints; ++q)
736 Blas::Dgemm(
'n',
'n', numPoints, numModes, numPoints, 1.0, D,
737 numPoints,
m_bdata.data(), numPoints, 0.0,
745 "Fourier modes should be a factor of 2");
747 for (i = 0; i < numPoints; ++i)
755 for (
p = 1;
p < numModes / 2; ++
p)
757 for (i = 0; i < numPoints; ++i)
759 m_bdata[2 *
p * numPoints + i] = cos(
p * M_PI * (z[i] + 1));
760 m_bdata[(2 *
p + 1) * numPoints + i] =
761 -sin(
p * M_PI * (z[i] + 1));
764 -
p * M_PI * sin(
p * M_PI * (z[i] + 1));
766 -
p * M_PI * cos(
p * M_PI * (z[i] + 1));
775 for (i = 0; i < numPoints; ++i)
777 m_bdata[i] = cos(M_PI * (z[i] + 1));
778 m_bdata[numPoints + i] = -sin(M_PI * (z[i] + 1));
780 m_dbdata[i] = -M_PI * sin(M_PI * (z[i] + 1));
781 m_dbdata[numPoints + i] = -M_PI * cos(M_PI * (z[i] + 1));
784 for (
p = 1;
p < numModes / 2; ++
p)
786 for (i = 0; i < numPoints; ++i)
788 m_bdata[2 *
p * numPoints + i] = 0.;
789 m_bdata[(2 *
p + 1) * numPoints + i] = 0.;
792 m_dbdata[(2 *
p + 1) * numPoints + i] = 0.;
800 m_dbdata[0] = -M_PI * sin(M_PI * z[0]);
805 m_bdata[0] = -sin(M_PI * z[0]);
806 m_dbdata[0] = -M_PI * cos(M_PI * z[0]);
813 for (
p = 0, scal = 1;
p < numModes; ++
p, mode += numPoints)
818 for (i = 0; i < numPoints; ++i)
823 scal *= 4 * (
p + 1) * (
p + 1) / (2 *
p + 2) / (2 *
p + 1);
827 Blas::Dgemm(
'n',
'n', numPoints, numModes, numPoints, 1.0, D,
828 numPoints,
m_bdata.data(), numPoints, 0.0,
835 int P = numModes - 1;
838 for (
int p = 0;
p <=
P; ++
p, mode += numPoints)
840 for (
int i = 0; i < numPoints; ++i)
842 mode[i] = pow(z[i],
p);
847 Blas::Dgemm(
'n',
'n', numPoints, numModes, numPoints, 1.0, D,
848 numPoints,
m_bdata.data(), numPoints, 0.0,
854 "not implemented at this time.");
863 bool returnval =
false;
871 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.
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.
static bool initBasisManager[]
Basis()
Private default constructor.
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.
Describes the specification for a Basis.
unsigned int m_nummodes
Expansion order.
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.
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 .
@ 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 .
The above copyright notice and this permission notice shall be included.
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