37#include <boost/core/ignore_unused.hpp> 
   72    const int nquad0 = 
m_base[0]->GetNumPoints();
 
   73    const int nquad1 = 
m_base[1]->GetNumPoints();
 
   74    const int nquad2 = 
m_base[2]->GetNumPoints();
 
   79    Vmath::Vcopy(nquad0 * nquad1 * nquad2, &inarray[0], 1, &wsp[0], 1);
 
   81    if (out_dx.size() > 0)
 
   85        Blas::Dgemm(
'N', 
'N', nquad0, nquad1 * nquad2, nquad0, 1.0, D0, nquad0,
 
   86                    &wsp[0], nquad0, 0.0, &out_dx[0], nquad0);
 
   89    if (out_dy.size() > 0)
 
   92        for (
int j = 0; j < nquad2; ++j)
 
   95                        &wsp[j * nquad0 * nquad1], nquad0, D1, nquad1, 0.0,
 
   96                        &out_dy[j * nquad0 * nquad1], nquad0);
 
  100    if (out_dz.size() > 0)
 
  104        Blas::Dgemm(
'N', 
'T', nquad0 * nquad1, nquad2, nquad2, 1.0, &wsp[0],
 
  105                    nquad0 * nquad1, D2, nquad2, 0.0, &out_dz[0],
 
  116    bool doCheckCollDir0, 
bool doCheckCollDir1, 
bool doCheckCollDir2)
 
  119                            doCheckCollDir0, doCheckCollDir1, doCheckCollDir2);
 
  128    bool doCheckCollDir0, 
bool doCheckCollDir1, 
bool doCheckCollDir2)
 
  131                                   doCheckCollDir0, doCheckCollDir1,
 
  137    ASSERTL1((dir == 0) || (dir == 1) || (dir == 2), 
"Invalid direction.");
 
  139    const int nq0 = 
m_base[0]->GetNumPoints();
 
  140    const int nq1 = 
m_base[1]->GetNumPoints();
 
  141    const int nq2 = 
m_base[2]->GetNumPoints();
 
  142    const int nq  = nq0 * nq1 * nq2;
 
  143    const int nm0 = 
m_base[0]->GetNumModes();
 
  144    const int nm1 = 
m_base[1]->GetNumModes();
 
  157            for (
int i = 0; i < nq; i++)
 
  162                    m_base[2]->GetBdata(), tmp2, tmp5, wsp, 
false, 
true, 
true);
 
  168                    (*mat)(j, i) = tmp5[j];
 
  173            for (
int i = 0; i < nq; i++)
 
  178                    m_base[2]->GetBdata(), tmp2, tmp5, wsp, 
true, 
false, 
true);
 
  184                    (*mat)(j, i) = tmp5[j];
 
  189            for (
int i = 0; i < nq; i++)
 
  194                    m_base[2]->GetDbdata(), tmp2, tmp5, wsp, 
true, 
true, 
false);
 
  199                    (*mat)(j, i) = tmp5[j];
 
  225    const int nq0 = 
m_base[0]->GetNumPoints();
 
  226    const int nq1 = 
m_base[1]->GetNumPoints();
 
  227    const int nq2 = 
m_base[2]->GetNumPoints();
 
  233    for (
int i = 0; i < nq1 * nq2; ++i, ptr += nq0)
 
  235        wsp1[i] = StdExpansion::BaryEvaluate<0>(eta[0], ptr);
 
  238    for (
int i = 0; i < nq2; ++i)
 
  240        wsp2[i] = StdExpansion::BaryEvaluate<1>(eta[1], &wsp1[i * nq1]);
 
  243    return StdExpansion::BaryEvaluate<2>(eta[2], &wsp2[0]);
 
  252    int Qx = 
m_base[0]->GetNumPoints();
 
  253    int Qy = 
m_base[1]->GetNumPoints();
 
  254    int Qz = 
m_base[2]->GetNumPoints();
 
  264    interpolatingNodes = &I[0]->GetPtr()[0];
 
  266    Blas::Dgemv(
'T', Qx, Qy * Qz, 1.0, &physvals[0], Qx, &interpolatingNodes[0],
 
  267                1, 0.0, &sumFactorization_qr[0], 1);
 
  270    interpolatingNodes = &I[1]->GetPtr()[0];
 
  272    Blas::Dgemv(
'T', Qy, Qz, 1.0, &sumFactorization_qr[0], Qy,
 
  273                &interpolatingNodes[0], 1, 0.0, &sumFactorization_r[0], 1);
 
  276    interpolatingNodes = &I[2]->GetPtr()[0];
 
  277    value = 
Blas::Ddot(Qz, interpolatingNodes, 1, &sumFactorization_r[0], 1);
 
  285    std::array<NekDouble, 3> &firstOrderDerivs)
 
  287    boost::ignore_unused(coord, inarray, firstOrderDerivs);
 
  316        if (!(
m_base[0]->Collocation() && 
m_base[1]->Collocation() &&
 
  317              m_base[2]->Collocation()))
 
  348        int nquad0  = 
m_base[0]->GetNumPoints();
 
  349        int nquad1  = 
m_base[1]->GetNumPoints();
 
  350        int nquad2  = 
m_base[2]->GetNumPoints();
 
  351        int nmodes0 = 
m_base[0]->GetNumModes();
 
  352        int nmodes1 = 
m_base[1]->GetNumModes();
 
  353        int nmodes2 = 
m_base[2]->GetNumModes();
 
  354        int wspsize = max(nquad0 * nmodes2 * (nmodes1 + nquad1),
 
  355                          nquad0 * nquad1 * (nquad2 + nmodes0) +
 
  356                              nmodes0 * nmodes1 * nquad2);
 
  367        if (!(
m_base[0]->Collocation() && 
m_base[1]->Collocation() &&
 
  368              m_base[2]->Collocation()))
 
  379                                         wsp2, 
true, 
true, 
true);
 
  419    boost::ignore_unused(i);
 
  428    boost::ignore_unused(tid, maparray, signarray, traceOrient);
 
  442    if (maparray.size() != map2.size())
 
  447    for (
int i = 0; i < map2.size(); ++i)
 
  449        maparray[i] = map1[map2[i]];
 
  455    const int numpoints, 
const int nummodes)
 
  457    boost::ignore_unused(facedir);
 
  459    switch (faceDirBasisType)
 
  511    const int numpoints, 
const int nummodes)
 
  513    switch (faceDirBasisType)
 
  541                        numpoints, LibUtilities::eGaussRadauMAlpha1Beta0);
 
  569                        numpoints, LibUtilities::eGaussRadauMAlpha1Beta0);
 
  599                        numpoints, LibUtilities::eGaussRadauMAlpha1Beta0);
 
#define WARNINGL2(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....
Describes the specification for a Basis.
Defines a specification for a set of points.
virtual ~StdExpansion3D() override
virtual int v_GetNedges(void) const
void BwdTrans_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
virtual void v_GetEdgeInteriorToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards)
void IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
virtual void v_BwdTrans_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)=0
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
virtual int v_GetEdgeNcoeffs(const int i) const
virtual 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.
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray) override
Integrates the specified function over the domain.
virtual void v_GenStdMatBwdDeriv(const int dir, DNekMatSharedPtr &mat) override
virtual void v_GetTraceToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient, int P, int Q) override
virtual void v_IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)=0
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d1, Array< OneD, NekDouble > &outarray_d2, Array< OneD, NekDouble > &outarray_d3)
Calculate the 3D derivative in the local tensor/collapsed coordinate at the physical points.
The base class for all shapes.
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
void GetElmtTraceToTraceMap(const unsigned int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
void LaplacianMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, 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 GetTraceCoeffMap(const unsigned int traceid, Array< OneD, unsigned int > &maparray)
virtual void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
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 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 = alpha A x plus beta y where A[m x n].
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...
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
@ eModified_B
Principle Modified Functions .
@ P
Monomial polynomials .
@ eOrtho_A
Principle Orthogonal Functions .
@ eModified_C
Principle Modified Functions .
@ eGLL_Lagrange
Lagrange for SEM basis .
@ eOrtho_C
Principle Orthogonal Functions .
@ eModifiedPyr_C
Principle Modified Functions.
@ eOrtho_B
Principle Orthogonal Functions .
@ eModified_A
Principle Modified Functions .
@ eOrthoPyr_C
Principle Orthogonal Functions .
static const NekDouble kNekZeroTol
LibUtilities::BasisKey EvaluateTriFaceBasisKey(const int facedir, const LibUtilities::BasisType faceDirBasisType, const int numpoints, const int nummodes)
LibUtilities::BasisKey EvaluateQuadFaceBasisKey(const int facedir, const LibUtilities::BasisType faceDirBasisType, const int numpoints, const int nummodes)
The above copyright notice and this permission notice shall be included.
std::shared_ptr< DNekMat > DNekMatSharedPtr
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)