51    if (lhsPointsKey < rhsPointsKey)
 
   55    if (lhsPointsKey != rhsPointsKey)
 
  103    std::shared_ptr<Basis> returnval(
new Basis(bkey));
 
  104    returnval->Initialize();
 
  112             "Cannot call Basis initialisation with zero or negative order");
 
  114                                    "zero or negative numbers of points");
 
  147    std::shared_ptr<NekMatrix<NekDouble>> ftB(
 
  246            for (
p = 0; 
p < numModes; ++
p, mode += numPoints)
 
  250                scal = 
sqrt(0.5 * (2.0 * 
p + 1.0));
 
  251                for (i = 0; i < numPoints; ++i)
 
  258            Blas::Dgemm(
'n', 
'n', numPoints, numModes, numPoints, 1.0, D,
 
  259                        numPoints, 
m_bdata.data(), numPoints, 0.0,
 
  280            for (
size_t p = 0; 
p < numModes; ++
p)
 
  282                for (
size_t q = 0; 
q < numModes - 
p; ++
q, mode += numPoints)
 
  286                    for (
size_t j = 0; j < numPoints; ++j)
 
  289                            sqrt(
p + 
q + 1.0) * pow(0.5 * (1.0 - 
z[j]), 
p);
 
  295            Blas::Dgemm(
'n', 
'n', numPoints, numModes * (numModes + 1) / 2,
 
  296                        numPoints, 1.0, D, numPoints, 
m_bdata.data(), numPoints,
 
  314            size_t P = numModes - 1, Q = numModes - 1, R = numModes - 1;
 
  317            for (
size_t p = 0; 
p <= 
P; ++
p)
 
  319                for (
size_t q = 0; 
q <= Q - 
p; ++
q)
 
  321                    for (
size_t r = 0; r <= R - 
p - 
q; ++r, mode += numPoints)
 
  324                                         2 * 
p + 2 * 
q + 2.0, 0.0);
 
  325                        for (
size_t k = 0; k < numPoints; ++k)
 
  328                            mode[k] *= pow(0.5 * (1.0 - 
z[k]), 
p + 
q);
 
  331                            mode[k] *= 
sqrt(r + 
p + 
q + 1.5);
 
  339                        numModes * (numModes + 1) * (numModes + 2) / 6,
 
  340                        numPoints, 1.0, D, numPoints, 
m_bdata.data(), numPoints,
 
  366            size_t P = numModes - 1, Q = numModes - 1, R = numModes - 1;
 
  369            for (
size_t p = 0; 
p <= 
P; ++
p)
 
  371                for (
size_t q = 0; 
q <= Q; ++
q)
 
  373                    for (
size_t r = 0; r <= R - std::max(
p, 
q);
 
  374                         ++r, mode += numPoints)
 
  383                        size_t pq = std::max(
p + 
q, 
size_t(0));
 
  387                        for (
size_t k = 0; k < numPoints; ++k)
 
  390                            mode[k] *= pow(0.5 * (1.0 - 
z[k]), pq);
 
  393                            mode[k] *= 
sqrt(r + pq + 1.5);
 
  401                        numModes * (numModes + 1) * (numModes + 2) / 6,
 
  402                        numPoints, 1.0, D, numPoints, 
m_bdata.data(), numPoints,
 
  414            for (i = 0; i < numPoints; ++i)
 
  417                m_bdata[numPoints + i] = 0.5 * (1 + 
z[i]);
 
  420            mode = 
m_bdata.data() + 2 * numPoints;
 
  422            for (
p = 2; 
p < numModes; ++
p, mode += numPoints)
 
  427                for (i = 0; i < numPoints; ++i)
 
  434            Blas::Dgemm(
'n', 
'n', numPoints, numModes, numPoints, 1.0, D,
 
  435                        numPoints, 
m_bdata.data(), numPoints, 0.0,
 
  461            for (i = 0; i < numPoints; ++i)
 
  463                m_bdata[0 * numPoints + i] = 0.5 * (1 - 
z[i]);
 
  464                m_bdata[1 * numPoints + i] = 0.5 * (1 + 
z[i]);
 
  467            mode = 
m_bdata.data() + 2 * numPoints;
 
  469            for (
q = 2; 
q < numModes; ++
q, mode += numPoints)
 
  474                for (i = 0; i < numPoints; ++i)
 
  481            for (i = 0; i < numPoints; ++i)
 
  483                mode[i] = 0.5 * (1 - 
z[i]);
 
  488            for (
q = 2; 
q < numModes; ++
q, mode += numPoints)
 
  493                for (i = 0; i < numPoints; ++i)
 
  501            one_p_z     = 
m_bdata.data() + numPoints;
 
  503            for (
p = 2; 
p < numModes; ++
p)
 
  505                for (i = 0; i < numPoints; ++i)
 
  507                    mode[i] = 
m_bdata[i] * one_m_z_pow[i];
 
  513                for (
q = 1; 
q < numModes - 
p; ++
q, mode += numPoints)
 
  518                    for (i = 0; i < numPoints; ++i)
 
  520                        mode[i] *= one_m_z_pow[i] * one_p_z[i];
 
  525            Blas::Dgemm(
'n', 
'n', numPoints, numModes * (numModes + 1) / 2,
 
  526                        numPoints, 1.0, D, numPoints, 
m_bdata.data(), numPoints,
 
  570            for (
p = 0; 
p < numModes; ++
p)
 
  572                N = numPoints * (numModes - 
p) * (numModes - 
p + 1) / 2;
 
  575                B_offset += numPoints * (numModes - 
p);
 
  581                        numModes * (numModes + 1) * (numModes + 2) / 6,
 
  582                        numPoints, 1.0, D, numPoints, 
m_bdata.data(), numPoints,
 
  620            N = numPoints * (numModes) * (numModes + 1) / 2;
 
  624            B_offset += numPoints * (numModes);
 
  626            N = numPoints * (numModes - 1);
 
  632            N = numPoints * (numModes - 1) * (numModes) / 2;
 
  637            B_offset += numPoints * (numModes - 1);
 
  641            mode = 
m_bdata.data() + offset;
 
  643            for (
p = 2; 
p < numModes; ++
p)
 
  646                N = numPoints * (numModes - 
p);
 
  653                one_p_z = 
m_bdata.data() + numPoints;
 
  655                for (
q = 2; 
q < numModes; ++
q)
 
  658                    for (i = 0; i < numPoints; ++i)
 
  671                    for (
size_t r = 1; r < numModes - std::max(
p, 
q); ++r)
 
  674                                         2 * 
p + 2 * 
q - 3, 1.0);
 
  676                        for (i = 0; i < numPoints; ++i)
 
  678                            mode[i] *= one_m_z_pow[i] * one_p_z[i];
 
  687                        numModes * (numModes + 1) * (2 * numModes + 1) / 6,
 
  688                        numPoints, 1.0, D, numPoints, 
m_bdata.data(), numPoints,
 
  696            std::shared_ptr<Points<NekDouble>> 
m_points =
 
  700            for (
p = 0; 
p < numModes; ++
p, mode += numPoints)
 
  702                for (
q = 0; 
q < numPoints; ++
q)
 
  710            Blas::Dgemm(
'n', 
'n', numPoints, numModes, numPoints, 1.0, D,
 
  711                        numPoints, 
m_bdata.data(), numPoints, 0.0,
 
  720            std::shared_ptr<Points<NekDouble>> 
m_points =
 
  724            for (
p = 0; 
p < numModes; ++
p, mode += numPoints)
 
  726                for (
q = 0; 
q < numPoints; ++
q)
 
  734            Blas::Dgemm(
'n', 
'n', numPoints, numModes, numPoints, 1.0, D,
 
  735                        numPoints, 
m_bdata.data(), numPoints, 0.0,
 
  744                     "Fourier modes should be a factor of 2");
 
  746            for (i = 0; i < numPoints; ++i)
 
  754            for (
p = 1; 
p < numModes / 2; ++
p)
 
  756                for (i = 0; i < numPoints; ++i)
 
  758                    m_bdata[2 * 
p * numPoints + i] = cos(
p * M_PI * (
z[i] + 1));
 
  759                    m_bdata[(2 * 
p + 1) * numPoints + i] =
 
  760                        -sin(
p * M_PI * (
z[i] + 1));
 
  763                        -
p * M_PI * sin(
p * M_PI * (
z[i] + 1));
 
  765                        -
p * M_PI * cos(
p * M_PI * (
z[i] + 1));
 
  774            for (i = 0; i < numPoints; ++i)
 
  776                m_bdata[i]             = cos(M_PI * (
z[i] + 1));
 
  777                m_bdata[numPoints + i] = -sin(M_PI * (
z[i] + 1));
 
  779                m_dbdata[i]             = -M_PI * sin(M_PI * (
z[i] + 1));
 
  780                m_dbdata[numPoints + i] = -M_PI * cos(M_PI * (
z[i] + 1));
 
  783            for (
p = 1; 
p < numModes / 2; ++
p)
 
  785                for (i = 0; i < numPoints; ++i)
 
  787                    m_bdata[2 * 
p * numPoints + i]       = 0.;
 
  788                    m_bdata[(2 * 
p + 1) * numPoints + i] = 0.;
 
  791                    m_dbdata[(2 * 
p + 1) * numPoints + i] = 0.;
 
  817            for (
p = 0, scal = 1; 
p < numModes; ++
p, mode += numPoints)
 
  822                for (i = 0; i < numPoints; ++i)
 
  827                scal *= 4 * (
p + 1) * (
p + 1) / (2 * 
p + 2) / (2 * 
p + 1);
 
  831            Blas::Dgemm(
'n', 
'n', numPoints, numModes, numPoints, 1.0, D,
 
  832                        numPoints, 
m_bdata.data(), numPoints, 0.0,
 
  841            for (
size_t p = 0; 
p < numModes; ++
p, mode += numPoints)
 
  843                for (
size_t i = 0; i < numPoints; ++i)
 
  845                    mode[i] = pow(
z[i], 
p);
 
  850            Blas::Dgemm(
'n', 
'n', numPoints, numModes, numPoints, 1.0, D,
 
  851                        numPoints, 
m_bdata.data(), numPoints, 0.0,
 
  858                                        "not implemented at this time.");
 
  867    bool returnval = 
false;
 
  875        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)
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