Nektar++
Public Member Functions | Protected Member Functions | Private Member Functions | List of all members
Nektar::StdRegions::StdPyrExp Class Reference

#include <StdPyrExp.h>

Inheritance diagram for Nektar::StdRegions::StdPyrExp:
[legend]

Public Member Functions

 StdPyrExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
 
 StdPyrExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc, NekDouble *coeffs, NekDouble *phys)
 
 StdPyrExp ()=default
 
 StdPyrExp (const StdPyrExp &T)=default
 
 ~StdPyrExp () override=default
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion3D
 StdExpansion3D (int numcoeffs, const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
 
 StdExpansion3D ()=default
 
 StdExpansion3D (const StdExpansion3D &T)=default
 
 ~StdExpansion3D () override=default
 
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. More...
 
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)
 
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)
 
int GetNedges () const
 return the number of edges in 3D expansion More...
 
int GetEdgeNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th edge. More...
 
void GetEdgeInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards)
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion
 StdExpansion ()
 Default Constructor. More...
 
 StdExpansion (const int numcoeffs, const int numbases, const LibUtilities::BasisKey &Ba=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bb=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bc=LibUtilities::NullBasisKey)
 Constructor. More...
 
 StdExpansion (const StdExpansion &T)
 Copy Constructor. More...
 
virtual ~StdExpansion ()
 Destructor. More...
 
int GetNumBases () const
 This function returns the number of 1D bases used in the expansion. More...
 
const Array< OneD, const LibUtilities::BasisSharedPtr > & GetBase () const
 This function gets the shared point to basis. More...
 
const LibUtilities::BasisSharedPtrGetBasis (int dir) const
 This function gets the shared point to basis in the dir direction. More...
 
int GetNcoeffs (void) const
 This function returns the total number of coefficients used in the expansion. More...
 
int GetTotPoints () const
 This function returns the total number of quadrature points used in the element. More...
 
LibUtilities::BasisType GetBasisType (const int dir) const
 This function returns the type of basis used in the dir direction. More...
 
int GetBasisNumModes (const int dir) const
 This function returns the number of expansion modes in the dir direction. More...
 
int EvalBasisNumModesMax (void) const
 This function returns the maximum number of expansion modes over all local directions. More...
 
LibUtilities::PointsType GetPointsType (const int dir) const
 This function returns the type of quadrature points used in the dir direction. More...
 
int GetNumPoints (const int dir) const
 This function returns the number of quadrature points in the dir direction. More...
 
const Array< OneD, const NekDouble > & GetPoints (const int dir) const
 This function returns a pointer to the array containing the quadrature points in dir direction. More...
 
int GetNverts () const
 This function returns the number of vertices of the expansion domain. More...
 
int GetTraceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th trace. More...
 
int GetTraceIntNcoeffs (const int i) const
 
int GetTraceNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th trace. More...
 
const LibUtilities::BasisKey GetTraceBasisKey (const int i, int k=-1) const
 This function returns the basis key belonging to the i-th trace. More...
 
LibUtilities::PointsKey GetTracePointsKey (const int i, int k=-1) const
 This function returns the basis key belonging to the i-th trace. More...
 
int NumBndryCoeffs (void) const
 
int NumDGBndryCoeffs (void) const
 
const LibUtilities::PointsKey GetNodalPointsKey () const
 This function returns the type of expansion Nodal point type if defined. More...
 
int GetNtraces () const
 Returns the number of trace elements connected to this element. More...
 
LibUtilities::ShapeType DetShapeType () const
 This function returns the shape of the expansion domain. More...
 
std::shared_ptr< StdExpansionGetStdExp () const
 
std::shared_ptr< StdExpansionGetLinStdExp (void) const
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion () const
 
bool IsNodalNonTensorialExp ()
 
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space. More...
 
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space. More...
 
void FwdTransBndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble Integral (const Array< OneD, const NekDouble > &inarray)
 This function integrates the specified function over the domain. More...
 
void FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 This function fills the array outarray with the mode-th mode of the expansion. More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 this function calculates the inner product of a given function f with the different modes of the expansion More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check)
 
void IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
int GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void SetElmtId (const int id)
 Set the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2=NullNekDouble1DArray, Array< OneD, NekDouble > &coords_3=NullNekDouble1DArray)
 this function returns the physical coordinates of the quadrature points of the expansion More...
 
void GetCoord (const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
 given the coordinates of a point of the element in the local collapsed coordinate system, this function calculates the physical coordinates of the point More...
 
DNekMatSharedPtr GetStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr GetStdStaticCondMatrix (const StdMatrixKey &mkey)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
void DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
int CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
NekDouble StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
int GetCoordim ()
 
void GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
void GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
int GetVertexMap (const int localVertexId, bool useCoeffPacking=false)
 
void GetTraceToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
 
void GetTraceCoeffMap (const unsigned int traceid, Array< OneD, unsigned int > &maparray)
 
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 GetTraceInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eForwards)
 
void GetTraceNumModes (const int tid, int &numModes0, int &numModes1, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2)
 
void MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
DNekMatSharedPtr CreateGeneralMatrix (const StdMatrixKey &mkey)
 this function generates the mass matrix \(\mathbf{M}[i][j] = \int \phi_i(\mathbf{x}) \phi_j(\mathbf{x}) d\mathbf{x}\) More...
 
void GeneralMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
void ExponentialFilter (Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
 
void LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
DNekMatSharedPtr GenMatrix (const StdMatrixKey &mkey)
 
void PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void PhysDeriv_s (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds)
 
void PhysDeriv_n (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn)
 
void PhysDirectionalDeriv (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &outarray)
 
void StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
 This function evaluates the first derivative of the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs, std::array< NekDouble, 6 > &secondOrderDerivs)
 
NekDouble PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluateBasis (const Array< OneD, const NekDouble > &coords, int mode)
 This function evaluates the basis function mode mode at a point coords of the domain. More...
 
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta. More...
 
void LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi)
 Convert local collapsed coordinates eta into local cartesian coordinate xi. More...
 
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual void v_DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
NekDouble Linf (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( L_\infty\) error \( |\epsilon|_\infty = \max |u - u_{exact}|\) where \( u_{exact}\) is given by the array sol. More...
 
NekDouble L2 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( L_2\) error, \( | \epsilon |_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 dx \right]^{1/2} d\xi_1 \) where \( u_{exact}\) is given by the array sol. More...
 
NekDouble H1 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( H^1\) error, \( | \epsilon |^1_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 + \nabla(u - u_{exact})\cdot\nabla(u - u_{exact})\cdot dx \right]^{1/2} d\xi_1 \) where \( u_{exact}\) is given by the array sol. More...
 
const LibUtilities::PointsKeyVector GetPointsKeys () const
 
DNekMatSharedPtr BuildInverseTransformationMatrix (const DNekScalMatSharedPtr &m_transformationmatrix)
 
void PhysInterpToSimplexEquiSpaced (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int npset=-1)
 This function performs an interpolation from the physical space points provided at input into an array of equispaced points which are not the collapsed coordinate. So for a tetrahedron you will only get a tetrahedral number of values. More...
 
void GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true)
 This function provides the connectivity of local simplices (triangles or tets) to connect the equispaced data points provided by PhysInterpToSimplexEquiSpaced. More...
 
void EquiSpacedToCoeffs (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs a projection/interpolation from the equispaced points sometimes used in post-processing onto the coefficient space. More...
 
template<class T >
std::shared_ptr< T > as ()
 
void IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
void GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 

Protected Member Functions

void v_PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2) override
 Calculate the derivative of the physical points. More...
 
void v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Calculate the derivative of the physical points in a given direction. More...
 
void v_StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2) override
 
void v_StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Backward transformation is evaluated at the quadrature points. More...
 
void v_BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
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) override
 
void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Forward transform from physical quadrature space stored in inarray and evaluate the expansion coefficients and store in outarray. More...
 
void v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Inner product of inarray over region with respect to the expansion basis m_base[0]->GetBdata(),m_base[1]->GetBdata(), m_base[2]->GetBdata() and return in outarray. More...
 
void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
 
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) override
 
void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta) override
 
void v_LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi) override
 
void v_GetCoords (Array< OneD, NekDouble > &xi_x, Array< OneD, NekDouble > &xi_y, Array< OneD, NekDouble > &xi_z) override
 
void v_FillMode (const int mode, Array< OneD, NekDouble > &outarray) override
 
void v_GetTraceNumModes (const int fid, int &numModes0, int &numModes1, Orientation faceOrient=eDir1FwdDir1_Dir2FwdDir2) override
 
NekDouble v_PhysEvaluateBasis (const Array< OneD, const NekDouble > &coords, int mode) final
 
NekDouble v_PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
 
int v_GetNverts () const override
 
int v_GetNedges () const override
 
int v_GetNtraces () const override
 
LibUtilities::ShapeType v_DetShapeType () const override
 
int v_NumBndryCoeffs () const override
 
int v_NumDGBndryCoeffs () const override
 
int v_GetTraceNcoeffs (const int i) const override
 
int v_GetTraceIntNcoeffs (const int i) const override
 
int v_GetTraceNumPoints (const int i) const override
 
int v_GetEdgeNcoeffs (const int i) const override
 
int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset) override
 
const LibUtilities::BasisKey v_GetTraceBasisKey (const int i, const int k) const override
 
int v_GetVertexMap (int localVertexId, bool useCoeffPacking=false) override
 
void v_GetInteriorMap (Array< OneD, unsigned int > &outarray) override
 
void v_GetBoundaryMap (Array< OneD, unsigned int > &outarray) override
 
void v_GetTraceCoeffMap (const unsigned int fid, Array< OneD, unsigned int > &maparray) override
 
void v_GetElmtTraceToTraceMap (const unsigned int fid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation faceOrient, int P, int Q) override
 
void v_GetEdgeInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2) override
 
void v_GetTraceInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2) override
 
DNekMatSharedPtr v_GenMatrix (const StdMatrixKey &mkey) override
 
DNekMatSharedPtr v_CreateStdMatrix (const StdMatrixKey &mkey) override
 
void v_SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey) override
 
void v_ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion3D
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. More...
 
NekDouble v_PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals) override
 
NekDouble v_PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
 
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_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 v_LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
void v_HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray) override
 Integrates the specified function over the domain. More...
 
virtual int v_GetNedges (void) const
 
virtual int v_GetEdgeNcoeffs (const int i) const
 
NekDouble BaryTensorDeriv (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
 
virtual void v_GetEdgeInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards)
 
void v_GetTraceToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient, int P, int Q) override
 
void v_GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat) override
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion
DNekMatSharedPtr CreateStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr CreateStdStaticCondMatrix (const StdMatrixKey &mkey)
 Create the static condensation of a matrix when using a boundary interior decomposition. More...
 
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase_SumFac (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GeneralMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
 
void LaplacianMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp_MatFree (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void HelmholtzMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
 
virtual NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
virtual void v_GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 
virtual void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals, NekDouble &deriv, NekDouble &deriv2)
 This function performs the barycentric interpolation of the polynomial stored in coord at a point physvals using barycentric interpolation weights in direction. More...
 
template<int DIR>
NekDouble BaryEvaluateBasis (const NekDouble &coord, const int &mode)
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals)
 Helper function to pass an unused value by reference into BaryEvaluate. More...
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals, NekDouble &deriv)
 

Private Member Functions

int GetMode (int I, int J, int K)
 Compute the mode number in the expansion for a particular tensorial combination. More...
 

Additional Inherited Members

- Protected Attributes inherited from Nektar::StdRegions::StdExpansion
Array< OneD, LibUtilities::BasisSharedPtrm_base
 
int m_elmt_id
 
int m_ncoeffs
 
LibUtilities::NekManager< StdMatrixKey, DNekMat, StdMatrixKey::opLessm_stdMatrixManager
 
LibUtilities::NekManager< StdMatrixKey, DNekBlkMat, StdMatrixKey::opLessm_stdStaticCondMatrixManager
 

Detailed Description

Definition at line 43 of file StdPyrExp.h.

Constructor & Destructor Documentation

◆ StdPyrExp() [1/4]

Nektar::StdRegions::StdPyrExp::StdPyrExp ( const LibUtilities::BasisKey Ba,
const LibUtilities::BasisKey Bb,
const LibUtilities::BasisKey Bc 
)

Definition at line 43 of file StdPyrExp.cpp.

47 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
48 3, Ba, Bb, Bc),
50 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
51 Ba, Bb, Bc)
52{
53
54 ASSERTL0(Ba.GetNumModes() <= Bc.GetNumModes(),
55 "order in 'a' direction is higher "
56 "than order in 'c' direction");
57 ASSERTL0(Bb.GetNumModes() <= Bc.GetNumModes(),
58 "order in 'b' direction is higher "
59 "than order in 'c' direction");
60 ASSERTL1(Bc.GetBasisType() == LibUtilities::eModifiedPyr_C ||
61 Bc.GetBasisType() == LibUtilities::eOrthoPyr_C,
62 "Expected basis type in 'c' direction to be ModifiedPyr_C or "
63 "OrthoPyr_C");
64}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
StdExpansion()
Default Constructor.
int getNumberOfCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:232
@ eModifiedPyr_C
Principle Modified Functions.
Definition: BasisType.h:53
@ eOrthoPyr_C
Principle Orthogonal Functions .
Definition: BasisType.h:51

References ASSERTL0, ASSERTL1, Nektar::LibUtilities::eModifiedPyr_C, Nektar::LibUtilities::eOrthoPyr_C, Nektar::LibUtilities::BasisKey::GetBasisType(), and Nektar::LibUtilities::BasisKey::GetNumModes().

◆ StdPyrExp() [2/4]

Nektar::StdRegions::StdPyrExp::StdPyrExp ( const LibUtilities::BasisKey Ba,
const LibUtilities::BasisKey Bb,
const LibUtilities::BasisKey Bc,
NekDouble coeffs,
NekDouble phys 
)

◆ StdPyrExp() [3/4]

Nektar::StdRegions::StdPyrExp::StdPyrExp ( )
default

◆ StdPyrExp() [4/4]

Nektar::StdRegions::StdPyrExp::StdPyrExp ( const StdPyrExp T)
default

◆ ~StdPyrExp()

Nektar::StdRegions::StdPyrExp::~StdPyrExp ( )
overridedefault

Member Function Documentation

◆ GetMode()

int Nektar::StdRegions::StdPyrExp::GetMode ( int  I,
int  J,
int  K 
)
private

Compute the mode number in the expansion for a particular tensorial combination.

Modes are numbered with the r index travelling fastest, followed by q and then p, and each q-r plane is of size

(R+1-p)*(Q+1) - l(l+1)/2 where l = max(0,Q-p)

For example, when P=2, Q=3 and R=4 the indexing inside each q-r plane (with r increasing upwards and q to the right) is:

p = 0: p = 1: p = 2:

4 3 8 17 21 2 7 11 16 20 24 29 32 35 1 6 10 13 15 19 23 26 28 31 34 37 0 5 9 12 14 18 22 25 27 30 33 36

Note that in this element, we must have that \( P,Q \leq R\).

Definition at line 1909 of file StdPyrExp.cpp.

1910{
1911 const int Q = m_base[1]->GetNumModes() - 1;
1912 const int R = m_base[2]->GetNumModes() - 1;
1913
1914 int i, l;
1915 int cnt = 0;
1916
1917 // Traverse to q-r plane number I
1918 for (i = 0; i < I; ++i)
1919 {
1920 // Size of triangle part
1921 l = max(0, Q - i);
1922
1923 // Size of rectangle part
1924 cnt += (R + 1 - i) * (Q + 1) - l * (l + 1) / 2;
1925 }
1926
1927 // Traverse to q column J (Pretend this is a face of width J)
1928 l = max(0, J - 1 - I);
1929 cnt += (R + 1 - I) * J - l * (l + 1) / 2;
1930
1931 // Traverse up stacks to K
1932 cnt += K;
1933
1934 return cnt;
1935}
Array< OneD, LibUtilities::BasisSharedPtr > m_base

References Nektar::StdRegions::StdExpansion::m_base.

Referenced by v_GetBoundaryMap(), v_GetEdgeInteriorToElementMap(), v_GetInteriorMap(), v_GetTraceCoeffMap(), v_GetTraceInteriorToElementMap(), and v_GetVertexMap().

◆ v_BwdTrans()

void Nektar::StdRegions::StdPyrExp::v_BwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Backward transformation is evaluated at the quadrature points.

\( u^{\delta} (\xi_{1i}, \xi_{2j}, \xi_{3k}) = \sum_{m(pqr)} \hat u_{pqr} \phi_{pqr} (\xi_{1i}, \xi_{2j}, \xi_{3k})\)

Backward transformation is three dimensional tensorial expansion

\( u (\xi_{1i}, \xi_{2j}, \xi_{3k}) = \sum_{p=0}^{Q_x} \psi_p^a (\xi_{1i}) \lbrace { \sum_{q=0}^{Q_y} \psi_{q}^a (\xi_{2j}) \lbrace { \sum_{r=0}^{Q_z} \hat u_{pqr} \psi_{pqr}^c (\xi_{3k}) \rbrace} \rbrace}. \)

And sumfactorizing step of the form is as:\ \ \( f_{pqr} (\xi_{3k}) = \sum_{r=0}^{Q_z} \hat u_{pqr} \psi_{pqr}^c (\xi_{3k}),\\ g_{p} (\xi_{2j}, \xi_{3k}) = \sum_{r=0}^{Q_y} \psi_{p}^a (\xi_{2j}) f_{pqr} (\xi_{3k}),\\ u(\xi_{1i}, \xi_{2j}, \xi_{3k}) = \sum_{p=0}^{Q_x} \psi_{p}^a (\xi_{1i}) g_{p} (\xi_{2j}, \xi_{3k}). \)

Implements Nektar::StdRegions::StdExpansion.

Definition at line 228 of file StdPyrExp.cpp.

230{
231 if (m_base[0]->Collocation() && m_base[1]->Collocation() &&
232 m_base[2]->Collocation())
233 {
235 m_base[2]->GetNumPoints(),
236 inarray, 1, outarray, 1);
237 }
238 else
239 {
240 StdPyrExp::v_BwdTrans_SumFac(inarray, outarray);
241 }
242}
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:218
void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdPyrExp.cpp:247
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References Nektar::StdRegions::StdExpansion::GetNumPoints(), Nektar::StdRegions::StdExpansion::m_base, v_BwdTrans_SumFac(), and Vmath::Vcopy().

Referenced by v_FillMode().

◆ v_BwdTrans_SumFac()

void Nektar::StdRegions::StdPyrExp::v_BwdTrans_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Sum-factorisation implementation of the BwdTrans operation.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 247 of file StdPyrExp.cpp.

249{
250 int nquad0 = m_base[0]->GetNumPoints();
251 int nquad1 = m_base[1]->GetNumPoints();
252 int nquad2 = m_base[2]->GetNumPoints();
253 int order0 = m_base[0]->GetNumModes();
254 int order1 = m_base[1]->GetNumModes();
255
256 Array<OneD, NekDouble> wsp(nquad2 * order0 * order1 +
257 nquad2 * nquad1 * nquad0);
258
259 v_BwdTrans_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetBdata(),
260 m_base[2]->GetBdata(), inarray, outarray, wsp, true,
261 true, true);
262}
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) override
Definition: StdPyrExp.cpp:264

References Nektar::StdRegions::StdExpansion::m_base, and v_BwdTrans_SumFacKernel().

Referenced by v_BwdTrans().

◆ v_BwdTrans_SumFacKernel()

void Nektar::StdRegions::StdPyrExp::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 
)
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion3D.

Definition at line 264 of file StdPyrExp.cpp.

273{
274 int nquad0 = m_base[0]->GetNumPoints();
275 int nquad1 = m_base[1]->GetNumPoints();
276 int nquad2 = m_base[2]->GetNumPoints();
277
278 int order0 = m_base[0]->GetNumModes();
279 int order1 = m_base[1]->GetNumModes();
280 int order2 = m_base[2]->GetNumModes();
281
282 Array<OneD, NekDouble> tmp = wsp;
283 Array<OneD, NekDouble> tmp1 = tmp + nquad2 * order0 * order1;
284
285 int i, j, mode, mode1, cnt;
286
287 // Perform summation over '2' direction
288 mode = mode1 = cnt = 0;
289 for (i = 0; i < order0; ++i)
290 {
291 for (j = 0; j < order1; ++j, ++cnt)
292 {
293 int ijmax = max(i, j);
294 Blas::Dgemv('N', nquad2, order2 - ijmax, 1.0,
295 base2.get() + mode * nquad2, nquad2,
296 inarray.get() + mode1, 1, 0.0, tmp.get() + cnt * nquad2,
297 1);
298 mode += order2 - ijmax;
299 mode1 += order2 - ijmax;
300 }
301 // increment mode in case order1!=order2
302 for (j = order1; j < order2; ++j)
303 {
304 int ijmax = max(i, j);
305 mode += order2 - ijmax;
306 }
307 }
308
309 // fix for modified basis by adding split of top singular
310 // vertex mode - currently (1+c)/2 x (1-b)/2 x (1-a)/2
311 // component is evaluated
313 {
314
315 // Not sure why we could not use basis as 1.0
316 // top singular vertex - (1+c)/2 x (1+b)/2 x (1-a)/2 component
317 Blas::Daxpy(nquad2, inarray[1], base2.get() + nquad2, 1,
318 &tmp[0] + nquad2, 1);
319
320 // top singular vertex - (1+c)/2 x (1-b)/2 x (1+a)/2 component
321 Blas::Daxpy(nquad2, inarray[1], base2.get() + nquad2, 1,
322 &tmp[0] + order1 * nquad2, 1);
323
324 // top singular vertex - (1+c)/2 x (1+b)/2 x (1+a)/2 component
325 Blas::Daxpy(nquad2, inarray[1], base2.get() + nquad2, 1,
326 &tmp[0] + order1 * nquad2 + nquad2, 1);
327 }
328
329 // Perform summation over '1' direction
330 mode = 0;
331 for (i = 0; i < order0; ++i)
332 {
333 Blas::Dgemm('N', 'T', nquad1, nquad2, order1, 1.0, base1.get(), nquad1,
334 tmp.get() + mode * nquad2, nquad2, 0.0,
335 tmp1.get() + i * nquad1 * nquad2, nquad1);
336 mode += order1;
337 }
338
339 // Perform summation over '0' direction
340 Blas::Dgemm('N', 'T', nquad0, nquad1 * nquad2, order0, 1.0, base0.get(),
341 nquad0, tmp1.get(), nquad1 * nquad2, 0.0, outarray.get(),
342 nquad0);
343}
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:156
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].
Definition: Blas.hpp:211
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...
Definition: Blas.hpp:383
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
Definition: Blas.hpp:135
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48

References Blas::Daxpy(), Blas::Dgemm(), Blas::Dgemv(), Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), and Nektar::StdRegions::StdExpansion::m_base.

Referenced by v_BwdTrans_SumFac().

◆ v_CalcNumberOfCoefficients()

int Nektar::StdRegions::StdPyrExp::v_CalcNumberOfCoefficients ( const std::vector< unsigned int > &  nummodes,
int &  modes_offset 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1028 of file StdPyrExp.cpp.

1030{
1032 nummodes[modes_offset], nummodes[modes_offset + 1],
1033 nummodes[modes_offset + 2]);
1034
1035 modes_offset += 3;
1036 return nmodes;
1037}

References Nektar::LibUtilities::StdPyrData::getNumberOfCoefficients().

◆ v_CreateStdMatrix()

DNekMatSharedPtr Nektar::StdRegions::StdPyrExp::v_CreateStdMatrix ( const StdMatrixKey mkey)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PyrExp.

Definition at line 1880 of file StdPyrExp.cpp.

1881{
1882 return v_GenMatrix(mkey);
1883}
DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override
Definition: StdPyrExp.cpp:1875

References v_GenMatrix().

◆ v_DetShapeType()

LibUtilities::ShapeType Nektar::StdRegions::StdPyrExp::v_DetShapeType ( ) const
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 871 of file StdPyrExp.cpp.

872{
874}

References Nektar::LibUtilities::ePyramid.

◆ v_FillMode()

void Nektar::StdRegions::StdPyrExp::v_FillMode ( const int  mode,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 752 of file StdPyrExp.cpp.

753{
754 Array<OneD, NekDouble> tmp(m_ncoeffs, 0.0);
755 tmp[mode] = 1.0;
756 v_BwdTrans(tmp, outarray);
757}
void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Backward transformation is evaluated at the quadrature points.
Definition: StdPyrExp.cpp:228

References Nektar::StdRegions::StdExpansion::m_ncoeffs, and v_BwdTrans().

◆ v_FwdTrans()

void Nektar::StdRegions::StdPyrExp::v_FwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Forward transform from physical quadrature space stored in inarray and evaluate the expansion coefficients and store in outarray.

Inputs:

  • inarray: array of physical quadrature points to be transformed

Outputs:

  • outarray: updated array of expansion coefficients.

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PyrExp.

Definition at line 358 of file StdPyrExp.cpp.

360{
361 v_IProductWRTBase(inarray, outarray);
362
363 // get Mass matrix inverse
364 StdMatrixKey imasskey(eInvMass, DetShapeType(), *this);
365 DNekMatSharedPtr imatsys = GetStdMatrix(imasskey);
366
367 // copy inarray in case inarray == outarray
368 DNekVec in(m_ncoeffs, outarray);
369 DNekVec out(m_ncoeffs, outarray, eWrapper);
370
371 out = (*imatsys) * in;
372}
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:603
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:367
void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Inner product of inarray over region with respect to the expansion basis m_base[0]->GetBdata(),...
Definition: StdPyrExp.cpp:393
NekVector< NekDouble > DNekVec
Definition: NekTypeDefs.hpp:48
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75

References Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::StdRegions::eInvMass, Nektar::eWrapper, Nektar::StdRegions::StdExpansion::GetStdMatrix(), Nektar::StdRegions::StdExpansion::m_ncoeffs, and v_IProductWRTBase().

◆ v_GenMatrix()

DNekMatSharedPtr Nektar::StdRegions::StdPyrExp::v_GenMatrix ( const StdMatrixKey mkey)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PyrExp.

Definition at line 1875 of file StdPyrExp.cpp.

1876{
1877 return CreateGeneralMatrix(mkey);
1878}
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix

References Nektar::StdRegions::StdExpansion::CreateGeneralMatrix().

Referenced by v_CreateStdMatrix().

◆ v_GetBoundaryMap()

void Nektar::StdRegions::StdPyrExp::v_GetBoundaryMap ( Array< OneD, unsigned int > &  outarray)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1137 of file StdPyrExp.cpp.

1138{
1141 "BasisType is not a boundary interior form");
1144 "BasisType is not a boundary interior form");
1147 "BasisType is not a boundary interior form");
1148
1149 int P = m_base[0]->GetNumModes() - 1, p;
1150 int Q = m_base[1]->GetNumModes() - 1, q;
1151 int R = m_base[2]->GetNumModes() - 1, r;
1152 int idx = 0;
1153
1154 int nBnd = NumBndryCoeffs();
1155
1156 if (maparray.size() != nBnd)
1157 {
1158 maparray = Array<OneD, unsigned int>(nBnd);
1159 }
1160
1161 // Loop over all boundary modes (in ascending order).
1162 for (p = 0; p <= P; ++p)
1163 {
1164 // First two q-r planes are entirely boundary modes.
1165 if (p <= 1)
1166 {
1167 for (q = 0; q <= Q; ++q)
1168 {
1169 int maxpq = max(p, q);
1170 for (r = 0; r <= R - maxpq; ++r)
1171 {
1172 maparray[idx++] = GetMode(p, q, r);
1173 }
1174 }
1175 }
1176 else
1177 {
1178 // Remaining q-r planes contain boundary modes on the two
1179 // front and back sides and edges 0 2.
1180 for (q = 0; q <= Q; ++q)
1181 {
1182 if (q <= 1)
1183 {
1184 for (r = 0; r <= R - p; ++r)
1185 {
1186 maparray[idx++] = GetMode(p, q, r);
1187 }
1188 }
1189 else
1190 {
1191 maparray[idx++] = GetMode(p, q, 0);
1192 }
1193 }
1194 }
1195 }
1196}
int GetMode(int I, int J, int K)
Compute the mode number in the expansion for a particular tensorial combination.
Definition: StdPyrExp.cpp:1909
@ P
Monomial polynomials .
Definition: BasisType.h:62
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:56
std::vector< double > q(NPUPPER *NPUPPER)

References ASSERTL1, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModifiedPyr_C, Nektar::StdRegions::StdExpansion::GetBasisType(), GetMode(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::NumBndryCoeffs(), Nektar::LibUtilities::P, CellMLToNektar.cellml_metadata::p, and Nektar::UnitTests::q().

◆ v_GetCoords()

void Nektar::StdRegions::StdPyrExp::v_GetCoords ( Array< OneD, NekDouble > &  xi_x,
Array< OneD, NekDouble > &  xi_y,
Array< OneD, NekDouble > &  xi_z 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PyrExp.

Definition at line 683 of file StdPyrExp.cpp.

686{
687 Array<OneD, const NekDouble> etaBar_x = m_base[0]->GetZ();
688 Array<OneD, const NekDouble> eta_y = m_base[1]->GetZ();
689 Array<OneD, const NekDouble> eta_z = m_base[2]->GetZ();
690 int Qx = GetNumPoints(0);
691 int Qy = GetNumPoints(1);
692 int Qz = GetNumPoints(2);
693
694 // Convert collapsed coordinates into cartesian coordinates: eta --> xi
695 for (int k = 0; k < Qz; ++k)
696 {
697 for (int j = 0; j < Qy; ++j)
698 {
699 for (int i = 0; i < Qx; ++i)
700 {
701 int s = i + Qx * (j + Qy * k);
702
703 xi_z[s] = eta_z[k];
704 xi_y[s] = (1.0 + eta_y[j]) * (1.0 - eta_z[k]) / 2.0 - 1.0;
705 xi_x[s] = (1.0 + etaBar_x[i]) * (1.0 - eta_z[k]) / 2.0 - 1.0;
706 }
707 }
708 }
709}

References Nektar::StdRegions::StdExpansion::GetNumPoints(), and Nektar::StdRegions::StdExpansion::m_base.

◆ v_GetEdgeInteriorToElementMap()

void Nektar::StdRegions::StdPyrExp::v_GetEdgeInteriorToElementMap ( const int  tid,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray,
const Orientation  traceOrient = eDir1FwdDir1_Dir2FwdDir2 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 1595 of file StdPyrExp.cpp.

1598{
1599 int i;
1600 bool signChange;
1601 const int P = m_base[0]->GetNumModes() - 1;
1602 const int Q = m_base[1]->GetNumModes() - 1;
1603 const int R = m_base[2]->GetNumModes() - 1;
1604 const int nEdgeIntCoeffs = v_GetEdgeNcoeffs(eid) - 2;
1605
1606 if (maparray.size() != nEdgeIntCoeffs)
1607 {
1608 maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
1609 }
1610
1611 if (signarray.size() != nEdgeIntCoeffs)
1612 {
1613 signarray = Array<OneD, int>(nEdgeIntCoeffs, 1);
1614 }
1615 else
1616 {
1617 fill(signarray.get(), signarray.get() + nEdgeIntCoeffs, 1);
1618 }
1619
1620 // If edge is oriented backwards, change sign of modes which have
1621 // degree 2n+1, n >= 1.
1622 signChange = edgeOrient == eBackwards;
1623
1624 switch (eid)
1625 {
1626 case 0:
1627 for (i = 2; i <= P; ++i)
1628 {
1629 maparray[i - 2] = GetMode(i, 0, 0);
1630 }
1631 break;
1632
1633 case 1:
1634 for (i = 2; i <= Q; ++i)
1635 {
1636 maparray[i - 2] = GetMode(1, i, 0);
1637 }
1638 break;
1639 case 2:
1640 for (i = 2; i <= P; ++i)
1641 {
1642 maparray[i - 2] = GetMode(i, 1, 0);
1643 }
1644 break;
1645
1646 case 3:
1647 for (i = 2; i <= Q; ++i)
1648 {
1649 maparray[i - 2] = GetMode(0, i, 0);
1650 }
1651 break;
1652 case 4:
1653 for (i = 2; i <= R; ++i)
1654 {
1655 maparray[i - 2] = GetMode(0, 0, i);
1656 }
1657 break;
1658
1659 case 5:
1660 for (i = 1; i <= R - 1; ++i)
1661 {
1662 maparray[i - 1] = GetMode(1, 0, i);
1663 }
1664 break;
1665 case 6:
1666 for (i = 1; i <= R - 1; ++i)
1667 {
1668 maparray[i - 1] = GetMode(1, 1, i);
1669 }
1670 break;
1671
1672 case 7:
1673 for (i = 1; i <= R - 1; ++i)
1674 {
1675 maparray[i - 1] = GetMode(0, 1, i);
1676 }
1677 break;
1678 default:
1679 ASSERTL0(false, "Edge not defined.");
1680 break;
1681 }
1682
1683 if (signChange)
1684 {
1685 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1686 {
1687 signarray[i] = -1;
1688 }
1689 }
1690}
int v_GetEdgeNcoeffs(const int i) const override
Definition: StdPyrExp.cpp:976

References ASSERTL0, Nektar::StdRegions::eBackwards, GetMode(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LibUtilities::P, and v_GetEdgeNcoeffs().

◆ v_GetEdgeNcoeffs()

int Nektar::StdRegions::StdPyrExp::v_GetEdgeNcoeffs ( const int  i) const
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 976 of file StdPyrExp.cpp.

977{
978 ASSERTL2(i >= 0 && i <= 7, "edge id is out of range");
979
980 if (i == 0 || i == 2)
981 {
982 return GetBasisNumModes(0);
983 }
984 else if (i == 1 || i == 3)
985 {
986 return GetBasisNumModes(1);
987 }
988 else
989 {
990 return GetBasisNumModes(2);
991 }
992}
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:265
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:169

References ASSERTL2, and Nektar::StdRegions::StdExpansion::GetBasisNumModes().

Referenced by v_GetEdgeInteriorToElementMap().

◆ v_GetElmtTraceToTraceMap()

void Nektar::StdRegions::StdPyrExp::v_GetElmtTraceToTraceMap ( const unsigned int  fid,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray,
Orientation  faceOrient,
int  P,
int  Q 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1315 of file StdPyrExp.cpp.

1319{
1321 "Method only implemented if BasisType is identical"
1322 "in x and y directions");
1325 "Method only implemented for Modified_A BasisType"
1326 "(x and y direction) and ModifiedPyr_C BasisType (z "
1327 "direction)");
1328
1329 int i, j, k, p, r, nFaceCoeffs;
1330 int nummodesA = 0, nummodesB = 0;
1331
1332 int order0 = m_base[0]->GetNumModes();
1333 int order1 = m_base[1]->GetNumModes();
1334 int order2 = m_base[2]->GetNumModes();
1335
1336 switch (fid)
1337 {
1338 case 0:
1339 nummodesA = order0;
1340 nummodesB = order1;
1341 break;
1342 case 1:
1343 case 3:
1344 nummodesA = order0;
1345 nummodesB = order2;
1346 break;
1347 case 2:
1348 case 4:
1349 nummodesA = order1;
1350 nummodesB = order2;
1351 break;
1352 default:
1353 ASSERTL0(false, "fid must be between 0 and 4");
1354 }
1355
1356 if (P == -1)
1357 {
1358 P = nummodesA;
1359 Q = nummodesB;
1360 nFaceCoeffs = GetTraceNcoeffs(fid);
1361 }
1362 else if (fid > 0)
1363 {
1364 nFaceCoeffs = P * (2 * Q - P + 1) / 2;
1365 }
1366 else
1367 {
1368 nFaceCoeffs = P * Q;
1369 }
1370
1371 // Allocate the map array and sign array; set sign array to ones (+)
1372 if (maparray.size() != nFaceCoeffs)
1373 {
1374 maparray = Array<OneD, unsigned int>(nFaceCoeffs);
1375 }
1376
1377 if (signarray.size() != nFaceCoeffs)
1378 {
1379 signarray = Array<OneD, int>(nFaceCoeffs, 1);
1380 }
1381 else
1382 {
1383 fill(signarray.get(), signarray.get() + nFaceCoeffs, 1);
1384 }
1385
1386 // triangular faces
1387 if (fid > 0)
1388 {
1389 // zero signmap and set maparray to zero if elemental
1390 // modes are not as large as face modesl
1391 int idx = 0;
1392 int cnt = 0;
1393 int minPA = min(nummodesA, P);
1394 int minQB = min(nummodesB, Q);
1395
1396 for (j = 0; j < minPA; ++j)
1397 {
1398 // set maparray
1399 for (k = 0; k < minQB - j; ++k, ++cnt)
1400 {
1401 maparray[idx++] = cnt;
1402 }
1403
1404 cnt += nummodesB - minQB;
1405
1406 for (k = nummodesB - j; k < Q - j; ++k)
1407 {
1408 signarray[idx] = 0.0;
1409 maparray[idx++] = maparray[0];
1410 }
1411 }
1412#if 0 // not required?
1413
1414 for (j = minPA; j < nummodesA; ++j)
1415 {
1416 // set maparray
1417 for (k = 0; k < minQB-j; ++k, ++cnt)
1418 {
1419 maparray[idx++] = cnt;
1420 }
1421
1422 cnt += nummodesB-minQB;
1423
1424 for (k = nummodesB-j; k < Q-j; ++k)
1425 {
1426 signarray[idx] = 0.0;
1427 maparray[idx++] = maparray[0];
1428 }
1429 }
1430#endif
1431 for (j = nummodesA; j < P; ++j)
1432 {
1433 for (k = 0; k < Q - j; ++k)
1434 {
1435 signarray[idx] = 0.0;
1436 maparray[idx++] = maparray[0];
1437 }
1438 }
1439
1440 // Triangles only have one possible orientation (base
1441 // direction reversed); swap edge modes.
1442 if (faceOrient == eDir1BwdDir1_Dir2FwdDir2)
1443 {
1444 swap(maparray[0], maparray[Q]);
1445 for (i = 1; i < Q - 1; ++i)
1446 {
1447 swap(maparray[i + 1], maparray[Q + i]);
1448 }
1449
1450 idx = 0;
1451 for (p = 0; p < P; ++p)
1452 {
1453 for (r = 0; r < Q - p; ++r, idx++)
1454 {
1455 if (p > 1)
1456 {
1457 signarray[idx] = p % 2 ? -1 : 1;
1458 }
1459 }
1460 }
1461 }
1462 }
1463 else
1464 {
1465
1466 // Set up an array indexing for quads, since the ordering may need
1467 // to be transposed.
1468 Array<OneD, int> arrayindx(nFaceCoeffs, -1);
1469
1470 for (i = 0; i < Q; i++)
1471 {
1472 for (j = 0; j < P; j++)
1473 {
1474 if (faceOrient < eDir1FwdDir2_Dir2FwdDir1)
1475 {
1476 arrayindx[i * P + j] = i * P + j;
1477 }
1478 else
1479 {
1480 arrayindx[i * P + j] = j * Q + i;
1481 }
1482 }
1483 }
1484
1485 // zero signmap and set maparray to zero if elemental
1486 // modes are not as large as face modesl
1487 for (j = 0; j < P; ++j)
1488 {
1489 // set up default maparray
1490 for (k = 0; k < Q; k++)
1491 {
1492 maparray[arrayindx[j + k * P]] = j + k * nummodesA;
1493 }
1494
1495 for (k = nummodesB; k < Q; ++k)
1496 {
1497 signarray[arrayindx[j + k * P]] = 0.0;
1498 maparray[arrayindx[j + k * P]] = maparray[0];
1499 }
1500 }
1501
1502 for (j = nummodesA; j < P; ++j)
1503 {
1504 for (k = 0; k < Q; ++k)
1505 {
1506 signarray[arrayindx[j + k * P]] = 0.0;
1507 maparray[arrayindx[j + k * P]] = maparray[0];
1508 }
1509 }
1510
1511 // The code below is exactly the same as that taken from
1512 // StdHexExp and reverses the 'b' and 'a' directions as
1513 // appropriate (1st and 2nd if statements respectively) in
1514 // quadrilateral faces.
1515 if (faceOrient == eDir1FwdDir1_Dir2BwdDir2 ||
1516 faceOrient == eDir1BwdDir1_Dir2BwdDir2 ||
1517 faceOrient == eDir1BwdDir2_Dir2FwdDir1 ||
1518 faceOrient == eDir1BwdDir2_Dir2BwdDir1)
1519 {
1520 if (faceOrient < eDir1FwdDir2_Dir2FwdDir1)
1521 {
1522 for (i = 3; i < Q; i += 2)
1523 {
1524 for (j = 0; j < P; j++)
1525 {
1526 signarray[arrayindx[i * P + j]] *= -1;
1527 }
1528 }
1529
1530 for (i = 0; i < P; i++)
1531 {
1532 swap(maparray[i], maparray[i + P]);
1533 swap(signarray[i], signarray[i + P]);
1534 }
1535 }
1536 else
1537 {
1538 for (i = 0; i < Q; i++)
1539 {
1540 for (j = 3; j < P; j += 2)
1541 {
1542 signarray[arrayindx[i * P + j]] *= -1;
1543 }
1544 }
1545
1546 for (i = 0; i < Q; i++)
1547 {
1548 swap(maparray[i], maparray[i + Q]);
1549 swap(signarray[i], signarray[i + Q]);
1550 }
1551 }
1552 }
1553
1554 if (faceOrient == eDir1BwdDir1_Dir2FwdDir2 ||
1555 faceOrient == eDir1BwdDir1_Dir2BwdDir2 ||
1556 faceOrient == eDir1FwdDir2_Dir2BwdDir1 ||
1557 faceOrient == eDir1BwdDir2_Dir2BwdDir1)
1558 {
1559 if (faceOrient < eDir1FwdDir2_Dir2FwdDir1)
1560 {
1561 for (i = 0; i < Q; i++)
1562 {
1563 for (j = 3; j < P; j += 2)
1564 {
1565 signarray[arrayindx[i * P + j]] *= -1;
1566 }
1567 }
1568
1569 for (i = 0; i < Q; i++)
1570 {
1571 swap(maparray[i * P], maparray[i * P + 1]);
1572 swap(signarray[i * P], signarray[i * P + 1]);
1573 }
1574 }
1575 else
1576 {
1577 for (i = 3; i < Q; i += 2)
1578 {
1579 for (j = 0; j < P; j++)
1580 {
1581 signarray[arrayindx[i * P + j]] *= -1;
1582 }
1583 }
1584
1585 for (i = 0; i < P; i++)
1586 {
1587 swap(maparray[i * Q], maparray[i * Q + 1]);
1588 swap(signarray[i * Q], signarray[i * Q + 1]);
1589 }
1590 }
1591 }
1592 }
1593}
int GetTraceNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th trace.
Definition: StdExpansion.h:261

References ASSERTL0, ASSERTL1, Nektar::StdRegions::eDir1BwdDir1_Dir2BwdDir2, Nektar::StdRegions::eDir1BwdDir1_Dir2FwdDir2, Nektar::StdRegions::eDir1BwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1BwdDir2_Dir2FwdDir1, Nektar::StdRegions::eDir1FwdDir1_Dir2BwdDir2, Nektar::StdRegions::eDir1FwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1FwdDir2_Dir2FwdDir1, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModifiedPyr_C, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::GetTraceNcoeffs(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LibUtilities::P, and CellMLToNektar.cellml_metadata::p.

◆ v_GetInteriorMap()

void Nektar::StdRegions::StdPyrExp::v_GetInteriorMap ( Array< OneD, unsigned int > &  outarray)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1098 of file StdPyrExp.cpp.

1099{
1102 "BasisType is not a boundary interior form");
1105 "BasisType is not a boundary interior form");
1108 "BasisType is not a boundary interior form");
1109
1110 int P = m_base[0]->GetNumModes() - 1, p;
1111 int Q = m_base[1]->GetNumModes() - 1, q;
1112 int R = m_base[2]->GetNumModes() - 1, r;
1113
1114 int nIntCoeffs = m_ncoeffs - NumBndryCoeffs();
1115
1116 if (outarray.size() != nIntCoeffs)
1117 {
1118 outarray = Array<OneD, unsigned int>(nIntCoeffs);
1119 }
1120
1121 int idx = 0;
1122
1123 // Loop over all interior modes.
1124 for (p = 2; p <= P; ++p)
1125 {
1126 for (q = 2; q <= Q; ++q)
1127 {
1128 int maxpq = max(p, q);
1129 for (r = 1; r <= R - maxpq; ++r)
1130 {
1131 outarray[idx++] = GetMode(p, q, r);
1132 }
1133 }
1134 }
1135}

References ASSERTL1, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModifiedPyr_C, Nektar::StdRegions::StdExpansion::GetBasisType(), GetMode(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::StdRegions::StdExpansion::NumBndryCoeffs(), Nektar::LibUtilities::P, CellMLToNektar.cellml_metadata::p, and Nektar::UnitTests::q().

◆ v_GetNedges()

int Nektar::StdRegions::StdPyrExp::v_GetNedges ( void  ) const
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 861 of file StdPyrExp.cpp.

862{
863 return 8;
864}

◆ v_GetNtraces()

int Nektar::StdRegions::StdPyrExp::v_GetNtraces ( ) const
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 866 of file StdPyrExp.cpp.

867{
868 return 5;
869}

◆ v_GetNverts()

int Nektar::StdRegions::StdPyrExp::v_GetNverts ( ) const
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 856 of file StdPyrExp.cpp.

857{
858 return 5;
859}

◆ v_GetTraceBasisKey()

const LibUtilities::BasisKey Nektar::StdRegions::StdPyrExp::v_GetTraceBasisKey ( const int  i,
const int  k 
) const
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 994 of file StdPyrExp.cpp.

996{
997 ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
998 ASSERTL2(k >= 0 && k <= 1, "basis key id is out of range");
999
1000 switch (i)
1001 {
1002 case 0:
1003 {
1005 m_base[k]->GetNumPoints(),
1006 m_base[k]->GetNumModes());
1007 }
1008 case 1:
1009 case 3:
1010 {
1011 return EvaluateTriFaceBasisKey(k, m_base[2 * k]->GetBasisType(),
1012 m_base[2 * k]->GetNumPoints(),
1013 m_base[2 * k]->GetNumModes());
1014 }
1015 case 2:
1016 case 4:
1017 {
1018 return EvaluateTriFaceBasisKey(k, m_base[k + 1]->GetBasisType(),
1019 m_base[k + 1]->GetNumPoints(),
1020 m_base[k + 1]->GetNumModes());
1021 }
1022 }
1023
1024 // Should never get here.
1026}
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.
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)

References ASSERTL2, Nektar::StdRegions::EvaluateQuadFaceBasisKey(), Nektar::StdRegions::EvaluateTriFaceBasisKey(), Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::GetNumPoints(), Nektar::StdRegions::StdExpansion::m_base, and Nektar::LibUtilities::NullBasisKey().

◆ v_GetTraceCoeffMap()

void Nektar::StdRegions::StdPyrExp::v_GetTraceCoeffMap ( const unsigned int  fid,
Array< OneD, unsigned int > &  maparray 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1198 of file StdPyrExp.cpp.

1200{
1202 "Method only implemented if BasisType is identical"
1203 "in x and y directions");
1206 "Method only implemented for Modified_A BasisType"
1207 "(x and y direction) and ModifiedPyr_C BasisType (z "
1208 "direction)");
1209
1210 int p, q, r, P = 0, Q = 0, idx = 0;
1211
1212 int order0 = m_base[0]->GetNumModes();
1213 int order1 = m_base[1]->GetNumModes();
1214 int order2 = m_base[2]->GetNumModes();
1215
1216 switch (fid)
1217 {
1218 case 0:
1219 P = order0;
1220 Q = order1;
1221 break;
1222 case 1:
1223 case 3:
1224 P = order0;
1225 Q = order2;
1226 break;
1227 case 2:
1228 case 4:
1229 P = order1;
1230 Q = order2;
1231 break;
1232 default:
1233 ASSERTL0(false, "fid must be between 0 and 4");
1234 }
1235
1236 if (maparray.size() != P * Q)
1237 {
1238 maparray = Array<OneD, unsigned int>(P * Q);
1239 }
1240
1241 // Set up ordering inside each 2D face. Also for triangular faces,
1242 // populate signarray.
1243 switch (fid)
1244 {
1245 case 0: // Bottom quad
1246
1247 for (q = 0; q < Q; ++q)
1248 {
1249 for (p = 0; p < P; ++p)
1250 {
1251 maparray[q * P + p] = GetMode(p, q, 0);
1252 }
1253 }
1254 break;
1255
1256 case 1: // Front triangle
1257 for (p = 0; p < P; ++p)
1258 {
1259 for (r = 0; r < Q - p; ++r)
1260 {
1261 maparray[idx++] = GetMode(p, 0, r);
1262 }
1263 }
1264 break;
1265
1266 case 2: // Right triangle
1267 maparray[idx++] = GetMode(1, 0, 0);
1268 maparray[idx++] = GetMode(0, 0, 1);
1269 for (r = 1; r < Q - 1; ++r)
1270 {
1271 maparray[idx++] = GetMode(1, 0, r);
1272 }
1273
1274 for (q = 1; q < P; ++q)
1275 {
1276 for (r = 0; r < Q - q; ++r)
1277 {
1278 maparray[idx++] = GetMode(1, q, r);
1279 }
1280 }
1281 break;
1282
1283 case 3: // Rear triangle
1284 maparray[idx++] = GetMode(0, 1, 0);
1285 maparray[idx++] = GetMode(0, 0, 1);
1286 for (r = 1; r < Q - 1; ++r)
1287 {
1288 maparray[idx++] = GetMode(0, 1, r);
1289 }
1290
1291 for (p = 1; p < P; ++p)
1292 {
1293 for (r = 0; r < Q - p; ++r)
1294 {
1295 maparray[idx++] = GetMode(p, 1, r);
1296 }
1297 }
1298 break;
1299
1300 case 4: // Left triangle
1301 for (q = 0; q < P; ++q)
1302 {
1303 for (r = 0; r < Q - q; ++r)
1304 {
1305 maparray[idx++] = GetMode(0, q, r);
1306 }
1307 }
1308 break;
1309
1310 default:
1311 ASSERTL0(false, "Face to element map unavailable.");
1312 }
1313}

References ASSERTL0, ASSERTL1, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModifiedPyr_C, Nektar::StdRegions::StdExpansion::GetBasisType(), GetMode(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LibUtilities::P, CellMLToNektar.cellml_metadata::p, and Nektar::UnitTests::q().

◆ v_GetTraceInteriorToElementMap()

void Nektar::StdRegions::StdPyrExp::v_GetTraceInteriorToElementMap ( const int  tid,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray,
const Orientation  traceOrient = eDir1FwdDir1_Dir2FwdDir2 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1692 of file StdPyrExp.cpp.

1695{
1696 const int P = m_base[0]->GetNumModes() - 1;
1697 const int Q = m_base[1]->GetNumModes() - 1;
1698 const int R = m_base[2]->GetNumModes() - 1;
1699 const int nFaceIntCoeffs = v_GetTraceIntNcoeffs(fid);
1700 int p, q, r, idx = 0;
1701 int nummodesA = 0;
1702 int nummodesB = 0;
1703 int i, j;
1704
1705 if (maparray.size() != nFaceIntCoeffs)
1706 {
1707 maparray = Array<OneD, unsigned int>(nFaceIntCoeffs);
1708 }
1709
1710 if (signarray.size() != nFaceIntCoeffs)
1711 {
1712 signarray = Array<OneD, int>(nFaceIntCoeffs, 1);
1713 }
1714 else
1715 {
1716 fill(signarray.get(), signarray.get() + nFaceIntCoeffs, 1);
1717 }
1718
1719 // Set up an array indexing for quad faces, since the ordering may
1720 // need to be transposed depending on orientation.
1721 Array<OneD, int> arrayindx(nFaceIntCoeffs);
1722 if (fid == 0)
1723 {
1724 nummodesA = P - 1;
1725 nummodesB = Q - 1;
1726
1727 for (i = 0; i < nummodesB; i++)
1728 {
1729 for (j = 0; j < nummodesA; j++)
1730 {
1731 if (faceOrient < 9)
1732 {
1733 arrayindx[i * nummodesA + j] = i * nummodesA + j;
1734 }
1735 else
1736 {
1737 arrayindx[i * nummodesA + j] = j * nummodesB + i;
1738 }
1739 }
1740 }
1741 }
1742
1743 switch (fid)
1744 {
1745 case 0: // Bottom quad
1746 for (q = 2; q <= Q; ++q)
1747 {
1748 for (p = 2; p <= P; ++p)
1749 {
1750 maparray[arrayindx[(q - 2) * nummodesA + (p - 2)]] =
1751 GetMode(p, q, 0);
1752 }
1753 }
1754 break;
1755 case 1: // Front triangle
1756 for (p = 2; p <= P; ++p)
1757 {
1758 for (r = 1; r <= R - p; ++r)
1759 {
1760 if ((int)faceOrient == 7)
1761 {
1762 signarray[idx] = p % 2 ? -1 : 1;
1763 }
1764 maparray[idx++] = GetMode(p, 0, r);
1765 }
1766 }
1767 break;
1768 case 2: // Right triangle
1769 for (q = 2; q <= Q; ++q)
1770 {
1771 for (r = 1; r <= R - q; ++r)
1772 {
1773 if ((int)faceOrient == 7)
1774 {
1775 signarray[idx] = q % 2 ? -1 : 1;
1776 }
1777 maparray[idx++] = GetMode(1, q, r);
1778 }
1779 }
1780 break;
1781
1782 case 3: // Rear triangle
1783 for (p = 2; p <= P; ++p)
1784 {
1785 for (r = 1; r <= R - p; ++r)
1786 {
1787 if ((int)faceOrient == 7)
1788 {
1789 signarray[idx] = p % 2 ? -1 : 1;
1790 }
1791 maparray[idx++] = GetMode(p, 1, r);
1792 }
1793 }
1794 break;
1795
1796 case 4: // Left triangle
1797 for (q = 2; q <= Q; ++q)
1798 {
1799 for (r = 1; r <= R - q; ++r)
1800 {
1801 if ((int)faceOrient == 7)
1802 {
1803 signarray[idx] = q % 2 ? -1 : 1;
1804 }
1805 maparray[idx++] = GetMode(0, q, r);
1806 }
1807 }
1808 break;
1809 default:
1810 ASSERTL0(false, "Face interior map not available.");
1811 }
1812
1813 // Triangular faces are processed in the above switch loop; for
1814 // remaining quad faces, set up orientation if necessary.
1815 if (fid > 0)
1816 {
1817 return;
1818 }
1819
1820 if (faceOrient == 6 || faceOrient == 8 || faceOrient == 11 ||
1821 faceOrient == 12)
1822 {
1823 if (faceOrient < 9)
1824 {
1825 for (i = 1; i < nummodesB; i += 2)
1826 {
1827 for (j = 0; j < nummodesA; j++)
1828 {
1829 signarray[arrayindx[i * nummodesA + j]] *= -1;
1830 }
1831 }
1832 }
1833 else
1834 {
1835 for (i = 0; i < nummodesB; i++)
1836 {
1837 for (j = 1; j < nummodesA; j += 2)
1838 {
1839 signarray[arrayindx[i * nummodesA + j]] *= -1;
1840 }
1841 }
1842 }
1843 }
1844
1845 if (faceOrient == 7 || faceOrient == 8 || faceOrient == 10 ||
1846 faceOrient == 12)
1847 {
1848 if (faceOrient < 9)
1849 {
1850 for (i = 0; i < nummodesB; i++)
1851 {
1852 for (j = 1; j < nummodesA; j += 2)
1853 {
1854 signarray[arrayindx[i * nummodesA + j]] *= -1;
1855 }
1856 }
1857 }
1858 else
1859 {
1860 for (i = 1; i < nummodesB; i += 2)
1861 {
1862 for (j = 0; j < nummodesA; j++)
1863 {
1864 signarray[arrayindx[i * nummodesA + j]] *= -1;
1865 }
1866 }
1867 }
1868 }
1869}
int v_GetTraceIntNcoeffs(const int i) const override
Definition: StdPyrExp.cpp:936

References ASSERTL0, GetMode(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LibUtilities::P, CellMLToNektar.cellml_metadata::p, Nektar::UnitTests::q(), and v_GetTraceIntNcoeffs().

◆ v_GetTraceIntNcoeffs()

int Nektar::StdRegions::StdPyrExp::v_GetTraceIntNcoeffs ( const int  i) const
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 936 of file StdPyrExp.cpp.

937{
938 ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
939
940 int P = m_base[0]->GetNumModes() - 1;
941 int Q = m_base[1]->GetNumModes() - 1;
942 int R = m_base[2]->GetNumModes() - 1;
943
944 if (i == 0)
945 {
946 return (P - 1) * (Q - 1);
947 }
948 else if (i == 1 || i == 3)
949 {
950 return (P - 1) * (2 * (R - 1) - (P - 1) - 1) / 2;
951 }
952 else
953 {
954 return (Q - 1) * (2 * (R - 1) - (Q - 1) - 1) / 2;
955 }
956}

References ASSERTL2, Nektar::StdRegions::StdExpansion::m_base, and Nektar::LibUtilities::P.

Referenced by v_GetTraceInteriorToElementMap().

◆ v_GetTraceNcoeffs()

int Nektar::StdRegions::StdPyrExp::v_GetTraceNcoeffs ( const int  i) const
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 916 of file StdPyrExp.cpp.

917{
918 ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
919
920 if (i == 0)
921 {
922 return GetBasisNumModes(0) * GetBasisNumModes(1);
923 }
924 else if (i == 1 || i == 3)
925 {
926 int P = GetBasisNumModes(0) - 1, Q = GetBasisNumModes(2) - 1;
927 return Q + 1 + (P * (1 + 2 * Q - P)) / 2;
928 }
929 else
930 {
931 int P = GetBasisNumModes(1) - 1, Q = GetBasisNumModes(2) - 1;
932 return Q + 1 + (P * (1 + 2 * Q - P)) / 2;
933 }
934}

References ASSERTL2, Nektar::StdRegions::StdExpansion::GetBasisNumModes(), and Nektar::LibUtilities::P.

◆ v_GetTraceNumModes()

void Nektar::StdRegions::StdPyrExp::v_GetTraceNumModes ( const int  fid,
int &  numModes0,
int &  numModes1,
Orientation  faceOrient = eDir1FwdDir1_Dir2FwdDir2 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 759 of file StdPyrExp.cpp.

761{
762 int nummodes[3] = {m_base[0]->GetNumModes(), m_base[1]->GetNumModes(),
763 m_base[2]->GetNumModes()};
764 switch (fid)
765 {
766 // quad
767 case 0:
768 {
769 numModes0 = nummodes[0];
770 numModes1 = nummodes[1];
771 }
772 break;
773 case 1:
774 case 3:
775 {
776 numModes0 = nummodes[0];
777 numModes1 = nummodes[2];
778 }
779 break;
780 case 2:
781 case 4:
782 {
783 numModes0 = nummodes[1];
784 numModes1 = nummodes[2];
785 }
786 break;
787 }
788
789 if (faceOrient >= 9)
790 {
791 std::swap(numModes0, numModes1);
792 }
793}

References Nektar::StdRegions::StdExpansion::m_base.

◆ v_GetTraceNumPoints()

int Nektar::StdRegions::StdPyrExp::v_GetTraceNumPoints ( const int  i) const
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 958 of file StdPyrExp.cpp.

959{
960 ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
961
962 if (i == 0)
963 {
964 return m_base[0]->GetNumPoints() * m_base[1]->GetNumPoints();
965 }
966 else if (i == 1 || i == 3)
967 {
968 return m_base[0]->GetNumPoints() * m_base[2]->GetNumPoints();
969 }
970 else
971 {
972 return m_base[1]->GetNumPoints() * m_base[2]->GetNumPoints();
973 }
974}

References ASSERTL2, and Nektar::StdRegions::StdExpansion::m_base.

◆ v_GetVertexMap()

int Nektar::StdRegions::StdPyrExp::v_GetVertexMap ( int  localVertexId,
bool  useCoeffPacking = false 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1039 of file StdPyrExp.cpp.

1040{
1044 "Mapping not defined for this type of basis");
1045
1046 int l = 0;
1047
1048 if (useCoeffPacking == true) // follow packing of coefficients i.e q,r,p
1049 {
1050 switch (vId)
1051 {
1052 case 0:
1053 l = GetMode(0, 0, 0);
1054 break;
1055 case 1:
1056 l = GetMode(0, 0, 1);
1057 break;
1058 case 2:
1059 l = GetMode(0, 1, 0);
1060 break;
1061 case 3:
1062 l = GetMode(1, 0, 0);
1063 break;
1064 case 4:
1065 l = GetMode(1, 1, 0);
1066 break;
1067 default:
1068 ASSERTL0(false, "local vertex id must be between 0 and 4");
1069 }
1070 }
1071 else
1072 {
1073 switch (vId)
1074 {
1075 case 0:
1076 l = GetMode(0, 0, 0);
1077 break;
1078 case 1:
1079 l = GetMode(1, 0, 0);
1080 break;
1081 case 2:
1082 l = GetMode(1, 1, 0);
1083 break;
1084 case 3:
1085 l = GetMode(0, 1, 0);
1086 break;
1087 case 4:
1088 l = GetMode(0, 0, 1);
1089 break;
1090 default:
1091 ASSERTL0(false, "local vertex id must be between 0 and 4");
1092 }
1093 }
1094
1095 return l;
1096}

References ASSERTL0, ASSERTL1, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModifiedPyr_C, Nektar::StdRegions::StdExpansion::GetBasisType(), and GetMode().

◆ v_IProductWRTBase()

void Nektar::StdRegions::StdPyrExp::v_IProductWRTBase ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Inner product of inarray over region with respect to the expansion basis m_base[0]->GetBdata(),m_base[1]->GetBdata(), m_base[2]->GetBdata() and return in outarray.

Wrapper call to StdPyrExp::IProductWRTBase

Input:

  • inarray: array of function evaluated at the physical collocation points

Output:

  • outarray: array of inner product with respect to each basis over region

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PyrExp.

Definition at line 393 of file StdPyrExp.cpp.

395{
396 if (m_base[0]->Collocation() && m_base[1]->Collocation() &&
397 m_base[2]->Collocation())
398 {
399 v_MultiplyByStdQuadratureMetric(inarray, outarray);
400 }
401 else
402 {
403 StdPyrExp::v_IProductWRTBase_SumFac(inarray, outarray);
404 }
405}
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
Definition: StdPyrExp.cpp:407
void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdPyrExp.cpp:1937

References Nektar::StdRegions::StdExpansion::m_base, v_IProductWRTBase_SumFac(), and v_MultiplyByStdQuadratureMetric().

Referenced by v_FwdTrans().

◆ v_IProductWRTBase_SumFac()

void Nektar::StdRegions::StdPyrExp::v_IProductWRTBase_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
bool  multiplybyweights = true 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PyrExp.

Definition at line 407 of file StdPyrExp.cpp.

410{
411
412 int nquad1 = m_base[1]->GetNumPoints();
413 int nquad2 = m_base[2]->GetNumPoints();
414 int order0 = m_base[0]->GetNumModes();
415 int order1 = m_base[1]->GetNumModes();
416
417 Array<OneD, NekDouble> wsp(order0 * nquad2 * (nquad1 + order1));
418
419 if (multiplybyweights)
420 {
421 Array<OneD, NekDouble> tmp(inarray.size());
422
424
426 m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
427 tmp, outarray, wsp, true, true, true);
428 }
429 else
430 {
432 m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
433 inarray, outarray, wsp, true, true, true);
434 }
435}
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) override
Definition: StdPyrExp.cpp:437

References Nektar::StdRegions::StdExpansion::m_base, v_IProductWRTBase_SumFacKernel(), and v_MultiplyByStdQuadratureMetric().

Referenced by v_IProductWRTBase().

◆ v_IProductWRTBase_SumFacKernel()

void Nektar::StdRegions::StdPyrExp::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 
)
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion3D.

Definition at line 437 of file StdPyrExp.cpp.

446{
447 int nquad0 = m_base[0]->GetNumPoints();
448 int nquad1 = m_base[1]->GetNumPoints();
449 int nquad2 = m_base[2]->GetNumPoints();
450
451 int order0 = m_base[0]->GetNumModes();
452 int order1 = m_base[1]->GetNumModes();
453 int order2 = m_base[2]->GetNumModes();
454
455 ASSERTL1(wsp.size() >= nquad1 * nquad2 * order0 + nquad2 * order0 * order1,
456 "Insufficient workspace size");
457
458 Array<OneD, NekDouble> tmp1 = wsp;
459 Array<OneD, NekDouble> tmp2 = wsp + nquad1 * nquad2 * order0;
460
461 int i, j, mode, mode1, cnt;
462
463 // Inner product with respect to the '0' direction
464 Blas::Dgemm('T', 'N', nquad1 * nquad2, order0, nquad0, 1.0, inarray.get(),
465 nquad0, base0.get(), nquad0, 0.0, tmp1.get(), nquad1 * nquad2);
466
467 // Inner product with respect to the '1' direction
468 for (mode = i = 0; i < order0; ++i)
469 {
470 Blas::Dgemm('T', 'N', nquad2, order1, nquad1, 1.0,
471 tmp1.get() + i * nquad1 * nquad2, nquad1, base1.get(),
472 nquad1, 0.0, tmp2.get() + mode * nquad2, nquad2);
473 mode += order1;
474 }
475
476 // Inner product with respect to the '2' direction
477 mode = mode1 = cnt = 0;
478 for (i = 0; i < order0; ++i)
479 {
480 for (j = 0; j < order1; ++j, ++cnt)
481 {
482 int ijmax = max(i, j);
483
484 Blas::Dgemv('T', nquad2, order2 - ijmax, 1.0,
485 base2.get() + mode * nquad2, nquad2,
486 tmp2.get() + cnt * nquad2, 1, 0.0,
487 outarray.get() + mode1, 1);
488 mode += order2 - ijmax;
489 mode1 += order2 - ijmax;
490 }
491
492 // increment mode in case order1!=order2
493 for (j = order1; j < order2; ++j)
494 {
495 int ijmax = max(i, j);
496 mode += order2 - ijmax;
497 }
498 }
499
500 // fix for modified basis for top singular vertex component
501 // Already have evaluated (1+c)/2 (1-b)/2 (1-a)/2
503 {
504 // add in (1+c)/2 (1+b)/2 (1-a)/2 component
505 outarray[1] +=
506 Blas::Ddot(nquad2, base2.get() + nquad2, 1, &tmp2[nquad2], 1);
507
508 // add in (1+c)/2 (1-b)/2 (1+a)/2 component
509 outarray[1] += Blas::Ddot(nquad2, base2.get() + nquad2, 1,
510 &tmp2[nquad2 * order1], 1);
511
512 // add in (1+c)/2 (1+b)/2 (1+a)/2 component
513 outarray[1] += Blas::Ddot(nquad2, base2.get() + nquad2, 1,
514 &tmp2[nquad2 * order1 + nquad2], 1);
515 }
516}
static double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = .
Definition: Blas.hpp:163

References ASSERTL1, Blas::Ddot(), Blas::Dgemm(), Blas::Dgemv(), Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), and Nektar::StdRegions::StdExpansion::m_base.

Referenced by v_IProductWRTBase_SumFac().

◆ v_IProductWRTDerivBase()

void Nektar::StdRegions::StdPyrExp::v_IProductWRTDerivBase ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PyrExp.

Definition at line 518 of file StdPyrExp.cpp.

521{
522 StdPyrExp::v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
523}
void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdPyrExp.cpp:531

References v_IProductWRTDerivBase_SumFac().

◆ v_IProductWRTDerivBase_SumFac()

void Nektar::StdRegions::StdPyrExp::v_IProductWRTDerivBase_SumFac ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual
Parameters
inarrayFunction evaluated at physical collocation points.
outarrayInner product with respect to each basis function over the element.

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PyrExp.

Definition at line 531 of file StdPyrExp.cpp.

534{
535 int i;
536 int nquad0 = m_base[0]->GetNumPoints();
537 int nquad1 = m_base[1]->GetNumPoints();
538 int nquad2 = m_base[2]->GetNumPoints();
539 int nqtot = nquad0 * nquad1 * nquad2;
540
541 Array<OneD, NekDouble> gfac0(nquad0);
542 Array<OneD, NekDouble> gfac1(nquad1);
543 Array<OneD, NekDouble> gfac2(nquad2);
544 Array<OneD, NekDouble> tmp0(nqtot);
545
546 int order0 = m_base[0]->GetNumModes();
547 int order1 = m_base[1]->GetNumModes();
548
549 Array<OneD, NekDouble> wsp(nquad1 * nquad2 * order0 +
550 nquad2 * order0 * order1);
551
552 const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
553 const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
554 const Array<OneD, const NekDouble> &z2 = m_base[2]->GetZ();
555
556 // set up geometric factor: (1+z0)/2
557 for (i = 0; i < nquad0; ++i)
558 {
559 gfac0[i] = 0.5 * (1 + z0[i]);
560 }
561
562 // set up geometric factor: (1+z1)/2
563 for (i = 0; i < nquad1; ++i)
564 {
565 gfac1[i] = 0.5 * (1 + z1[i]);
566 }
567
568 // Set up geometric factor: 2/(1-z2)
569 for (i = 0; i < nquad2; ++i)
570 {
571 gfac2[i] = 2.0 / (1 - z2[i]);
572 }
573
574 // Derivative in first/second direction is always scaled as follows
575 const int nq01 = nquad0 * nquad1;
576 for (i = 0; i < nquad2; ++i)
577 {
578 Vmath::Smul(nq01, gfac2[i], &inarray[0] + i * nq01, 1,
579 &tmp0[0] + i * nq01, 1);
580 }
581
583
584 switch (dir)
585 {
586 case 0:
587 {
589 m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
590 m_base[2]->GetBdata(), tmp0, outarray, wsp, false, true, true);
591 break;
592 }
593 case 1:
594 {
596 m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
597 m_base[2]->GetBdata(), tmp0, outarray, wsp, true, false, true);
598 break;
599 }
600 case 2:
601 {
602 Array<OneD, NekDouble> tmp3(m_ncoeffs);
603 Array<OneD, NekDouble> tmp4(m_ncoeffs);
604
605 // Scale eta_1 derivative by gfac0
606 for (i = 0; i < nquad1 * nquad2; ++i)
607 {
608 Vmath::Vmul(nquad0, tmp0.get() + i * nquad0, 1, gfac0.get(), 1,
609 tmp0.get() + i * nquad0, 1);
610 }
612 m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
613 m_base[2]->GetBdata(), tmp0, tmp3, wsp, false, true, true);
614
615 // Scale eta_2 derivative by gfac1*gfac2
616 for (i = 0; i < nquad2; ++i)
617 {
618 Vmath::Smul(nq01, gfac2[i], &inarray[0] + i * nq01, 1,
619 &tmp0[0] + i * nq01, 1);
620 }
621 for (i = 0; i < nquad1 * nquad2; ++i)
622 {
623 Vmath::Smul(nquad0, gfac1[i % nquad1], &tmp0[0] + i * nquad0, 1,
624 &tmp0[0] + i * nquad0, 1);
625 }
626
629 m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
630 m_base[2]->GetBdata(), tmp0, tmp4, wsp, true, false, true);
631
632 v_MultiplyByStdQuadratureMetric(inarray, tmp0);
634 m_base[0]->GetBdata(), m_base[1]->GetBdata(),
635 m_base[2]->GetDbdata(), tmp0, outarray, wsp, true, true, false);
636
637 Vmath::Vadd(m_ncoeffs, &tmp3[0], 1, &outarray[0], 1, &outarray[0],
638 1);
639 Vmath::Vadd(m_ncoeffs, &tmp4[0], 1, &outarray[0], 1, &outarray[0],
640 1);
641 break;
642 }
643 default:
644 {
645 ASSERTL1(false, "input dir is out of range");
646 break;
647 }
648 }
649}
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)
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.
Definition: Vmath.hpp:72
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100

References ASSERTL1, Nektar::StdRegions::StdExpansion3D::IProductWRTBase_SumFacKernel(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Vmath::Smul(), v_MultiplyByStdQuadratureMetric(), Vmath::Vadd(), and Vmath::Vmul().

Referenced by v_IProductWRTDerivBase().

◆ v_LocCollapsedToLocCoord()

void Nektar::StdRegions::StdPyrExp::v_LocCollapsedToLocCoord ( const Array< OneD, const NekDouble > &  eta,
Array< OneD, NekDouble > &  xi 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 675 of file StdPyrExp.cpp.

677{
678 xi[0] = (1.0 + eta[0]) * (1.0 - eta[2]) * 0.5 - 1.0;
679 xi[1] = (1.0 + eta[1]) * (1.0 - eta[2]) * 0.5 - 1.0;
680 xi[2] = eta[2];
681}

◆ v_LocCoordToLocCollapsed()

void Nektar::StdRegions::StdPyrExp::v_LocCoordToLocCollapsed ( const Array< OneD, const NekDouble > &  xi,
Array< OneD, NekDouble > &  eta 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 655 of file StdPyrExp.cpp.

657{
658 NekDouble d2 = 1.0 - xi[2];
659 if (fabs(d2) < NekConstants::kNekZeroTol)
660 {
661 if (d2 >= 0.)
662 {
664 }
665 else
666 {
668 }
669 }
670 eta[2] = xi[2]; // eta_z = xi_z
671 eta[1] = 2.0 * (1.0 + xi[1]) / d2 - 1.0;
672 eta[0] = 2.0 * (1.0 + xi[0]) / d2 - 1.0;
673}
static const NekDouble kNekZeroTol
double NekDouble

References Nektar::NekConstants::kNekZeroTol.

◆ v_MultiplyByStdQuadratureMetric()

void Nektar::StdRegions::StdPyrExp::v_MultiplyByStdQuadratureMetric ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1937 of file StdPyrExp.cpp.

1940{
1941 int i, j;
1942
1943 int nquad0 = m_base[0]->GetNumPoints();
1944 int nquad1 = m_base[1]->GetNumPoints();
1945 int nquad2 = m_base[2]->GetNumPoints();
1946
1947 const Array<OneD, const NekDouble> &w0 = m_base[0]->GetW();
1948 const Array<OneD, const NekDouble> &w1 = m_base[1]->GetW();
1949 const Array<OneD, const NekDouble> &w2 = m_base[2]->GetW();
1950
1951 const Array<OneD, const NekDouble> &z2 = m_base[2]->GetZ();
1952
1953 // Multiply by integration constants in x-direction
1954 for (i = 0; i < nquad1 * nquad2; ++i)
1955 {
1956 Vmath::Vmul(nquad0, inarray.get() + i * nquad0, 1, w0.get(), 1,
1957 outarray.get() + i * nquad0, 1);
1958 }
1959
1960 // Multiply by integration constants in y-direction
1961 for (j = 0; j < nquad2; ++j)
1962 {
1963 for (i = 0; i < nquad1; ++i)
1964 {
1965 Blas::Dscal(nquad0, w1[i],
1966 &outarray[0] + i * nquad0 + j * nquad0 * nquad1, 1);
1967 }
1968 }
1969
1970 // Multiply by integration constants in z-direction; need to
1971 // incorporate factor [(1-eta_3)/2]^2 into weights, but only if
1972 // using GLL quadrature points.
1973 switch (m_base[2]->GetPointsType())
1974 {
1975 // (2,0) Jacobi inner product.
1976 case LibUtilities::eGaussRadauMAlpha2Beta0:
1977 for (i = 0; i < nquad2; ++i)
1978 {
1979 Blas::Dscal(nquad0 * nquad1, 0.25 * w2[i],
1980 &outarray[0] + i * nquad0 * nquad1, 1);
1981 }
1982 break;
1983
1984 default:
1985 for (i = 0; i < nquad2; ++i)
1986 {
1987 Blas::Dscal(nquad0 * nquad1,
1988 0.25 * (1 - z2[i]) * (1 - z2[i]) * w2[i],
1989 &outarray[0] + i * nquad0 * nquad1, 1);
1990 }
1991 break;
1992 }
1993}
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:205
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition: Blas.hpp:149

References Blas::Dscal(), Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::StdRegions::StdExpansion::m_base, and Vmath::Vmul().

Referenced by v_IProductWRTBase(), v_IProductWRTBase_SumFac(), and v_IProductWRTDerivBase_SumFac().

◆ v_NumBndryCoeffs()

int Nektar::StdRegions::StdPyrExp::v_NumBndryCoeffs ( ) const
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 876 of file StdPyrExp.cpp.

877{
880 "BasisType is not a boundary interior form");
883 "BasisType is not a boundary interior form");
886 "BasisType is not a boundary interior form");
887
888 int P = m_base[0]->GetNumModes();
889 int Q = m_base[1]->GetNumModes();
890 int R = m_base[2]->GetNumModes();
891
893}
int getNumberOfBndCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:259

References ASSERTL1, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModifiedPyr_C, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::LibUtilities::StdPyrData::getNumberOfBndCoefficients(), Nektar::StdRegions::StdExpansion::m_base, and Nektar::LibUtilities::P.

◆ v_NumDGBndryCoeffs()

int Nektar::StdRegions::StdPyrExp::v_NumDGBndryCoeffs ( ) const
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 895 of file StdPyrExp.cpp.

896{
899 "BasisType is not a boundary interior form");
902 "BasisType is not a boundary interior form");
905 "BasisType is not a boundary interior form");
906
907 int P = m_base[0]->GetNumModes() - 1;
908 int Q = m_base[1]->GetNumModes() - 1;
909 int R = m_base[2]->GetNumModes() - 1;
910
911 return (P + 1) * (Q + 1) // 1 rect. face on base
912 + 2 * (R + 1) + P * (1 + 2 * R - P) // 2 tri. (P,R) faces
913 + 2 * (R + 1) + Q * (1 + 2 * R - Q); // 2 tri. (Q,R) faces
914}

References ASSERTL1, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModifiedPyr_C, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::m_base, and Nektar::LibUtilities::P.

◆ v_PhysDeriv() [1/2]

void Nektar::StdRegions::StdPyrExp::v_PhysDeriv ( const Array< OneD, const NekDouble > &  u_physical,
Array< OneD, NekDouble > &  out_dxi1,
Array< OneD, NekDouble > &  out_dxi2,
Array< OneD, NekDouble > &  out_dxi3 
)
overrideprotectedvirtual

Calculate the derivative of the physical points.

The derivative is evaluated at the nodal physical points. Derivatives with respect to the local Cartesian coordinates.

\(\begin{Bmatrix} \frac {\partial} {\partial \xi_1} \\ \frac {\partial} {\partial \xi_2} \\ \frac {\partial} {\partial \xi_3} \end{Bmatrix} = \begin{Bmatrix} \frac 2 {(1-\eta_3)} \frac \partial {\partial \bar \eta_1} \\ \frac {\partial} {\partial \xi_2} \ \ \frac {(1 + \bar \eta_1)} {(1 - \eta_3)} \frac \partial {\partial \bar \eta_1} + \frac {\partial} {\partial \eta_3} \end{Bmatrix}\)

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PyrExp.

Definition at line 83 of file StdPyrExp.cpp.

87{
88 // PhysDerivative implementation based on Spen's book page 152.
89 int Qx = m_base[0]->GetNumPoints();
90 int Qy = m_base[1]->GetNumPoints();
91 int Qz = m_base[2]->GetNumPoints();
92
93 Array<OneD, NekDouble> dEta_bar1(Qx * Qy * Qz, 0.0);
94 Array<OneD, NekDouble> dXi2(Qx * Qy * Qz, 0.0);
95 Array<OneD, NekDouble> dEta3(Qx * Qy * Qz, 0.0);
96 PhysTensorDeriv(u_physical, dEta_bar1, dXi2, dEta3);
97
98 Array<OneD, const NekDouble> eta_x, eta_y, eta_z;
99 eta_x = m_base[0]->GetZ();
100 eta_y = m_base[1]->GetZ();
101 eta_z = m_base[2]->GetZ();
102
103 int i, j, k, n;
104
105 if (out_dxi1.size() > 0)
106 {
107 for (k = 0, n = 0; k < Qz; ++k)
108 {
109 NekDouble fac = 2.0 / (1.0 - eta_z[k]);
110 for (j = 0; j < Qy; ++j)
111 {
112 for (i = 0; i < Qx; ++i, ++n)
113 {
114 out_dxi1[n] = fac * dEta_bar1[n];
115 }
116 }
117 }
118 }
119
120 if (out_dxi2.size() > 0)
121 {
122 for (k = 0, n = 0; k < Qz; ++k)
123 {
124 NekDouble fac = 2.0 / (1.0 - eta_z[k]);
125 for (j = 0; j < Qy; ++j)
126 {
127 for (i = 0; i < Qx; ++i, ++n)
128 {
129 out_dxi2[n] = fac * dXi2[n];
130 }
131 }
132 }
133 }
134
135 if (out_dxi3.size() > 0)
136 {
137 for (k = 0, n = 0; k < Qz; ++k)
138 {
139 NekDouble fac = 1.0 / (1.0 - eta_z[k]);
140 for (j = 0; j < Qy; ++j)
141 {
142 NekDouble fac1 = (1.0 + eta_y[j]);
143 for (i = 0; i < Qx; ++i, ++n)
144 {
145 out_dxi3[n] = (1.0 + eta_x[i]) * fac * dEta_bar1[n] +
146 fac1 * fac * dXi2[n] + dEta3[n];
147 }
148 }
149 }
150 }
151}
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.

References Nektar::StdRegions::StdExpansion::m_base, and Nektar::StdRegions::StdExpansion3D::PhysTensorDeriv().

Referenced by v_PhysDeriv(), and v_StdPhysDeriv().

◆ v_PhysDeriv() [2/2]

void Nektar::StdRegions::StdPyrExp::v_PhysDeriv ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d0 
)
overrideprotectedvirtual

Calculate the derivative of the physical points in a given direction.

See also
StdRegions::StdExpansion::PhysDeriv

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 153 of file StdPyrExp.cpp.

156{
157 switch (dir)
158 {
159 case 0:
160 {
161 v_PhysDeriv(inarray, outarray, NullNekDouble1DArray,
163 break;
164 }
165
166 case 1:
167 {
168 v_PhysDeriv(inarray, NullNekDouble1DArray, outarray,
170 break;
171 }
172
173 case 2:
174 {
176 outarray);
177 break;
178 }
179
180 default:
181 {
182 ASSERTL1(false, "input dir is out of range");
183 }
184 break;
185 }
186}
void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2) override
Calculate the derivative of the physical points.
Definition: StdPyrExp.cpp:83
static Array< OneD, NekDouble > NullNekDouble1DArray

References ASSERTL1, Nektar::NullNekDouble1DArray, and v_PhysDeriv().

◆ v_PhysEvaluate()

NekDouble Nektar::StdRegions::StdPyrExp::v_PhysEvaluate ( const Array< OneD, NekDouble > &  coord,
const Array< OneD, const NekDouble > &  inarray,
std::array< NekDouble, 3 > &  firstOrderDerivs 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Reimplemented in Nektar::LocalRegions::PyrExp.

Definition at line 711 of file StdPyrExp.cpp.

714{
715 // Collapse coordinates
716 Array<OneD, NekDouble> coll(3, 0.0);
717 LocCoordToLocCollapsed(coord, coll);
718
719 // If near singularity do the old interpolation matrix method
720 if ((1 - coll[2]) < 1e-5)
721 {
722 int totPoints = GetTotPoints();
723 Array<OneD, NekDouble> EphysDeriv0(totPoints), EphysDeriv1(totPoints),
724 EphysDeriv2(totPoints);
725 PhysDeriv(inarray, EphysDeriv0, EphysDeriv1, EphysDeriv2);
726
727 Array<OneD, DNekMatSharedPtr> I(3);
728 I[0] = GetBase()[0]->GetI(coll);
729 I[1] = GetBase()[1]->GetI(coll + 1);
730 I[2] = GetBase()[2]->GetI(coll + 2);
731
732 firstOrderDerivs[0] = PhysEvaluate(I, EphysDeriv0);
733 firstOrderDerivs[1] = PhysEvaluate(I, EphysDeriv1);
734 firstOrderDerivs[2] = PhysEvaluate(I, EphysDeriv2);
735 return PhysEvaluate(I, inarray);
736 }
737
738 std::array<NekDouble, 3> interDeriv;
739 NekDouble val = StdExpansion3D::BaryTensorDeriv(coll, inarray, interDeriv);
740
741 NekDouble fac = 2.0 / (1.0 - coll[2]);
742
743 firstOrderDerivs[0] = fac * interDeriv[0];
744 firstOrderDerivs[1] = fac * interDeriv[1];
745 firstOrderDerivs[2] = ((1.0 + coll[0]) / (1.0 - coll[2])) * interDeriv[0] +
746 ((1.0 + coll[1]) / (1.0 - coll[2])) * interDeriv[1] +
747 interDeriv[2];
748
749 return val;
750}
NekDouble BaryTensorDeriv(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:134
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
const Array< OneD, const LibUtilities::BasisSharedPtr > & GetBase() const
This function gets the shared point to basis.
Definition: StdExpansion.h:100
NekDouble PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
This function evaluates the expansion at a single (arbitrary) point of the domain.
Definition: StdExpansion.h:919
void PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Definition: StdExpansion.h:849

References Nektar::StdRegions::StdExpansion3D::BaryTensorDeriv(), Nektar::StdRegions::StdExpansion::GetBase(), Nektar::StdRegions::StdExpansion::GetTotPoints(), Nektar::StdRegions::StdExpansion::LocCoordToLocCollapsed(), Nektar::StdRegions::StdExpansion::PhysDeriv(), and Nektar::StdRegions::StdExpansion::PhysEvaluate().

◆ v_PhysEvaluateBasis()

NekDouble Nektar::StdRegions::StdPyrExp::v_PhysEvaluateBasis ( const Array< OneD, const NekDouble > &  coords,
int  mode 
)
finalprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 795 of file StdPyrExp.cpp.

797{
798 Array<OneD, NekDouble> coll(3);
799 LocCoordToLocCollapsed(coords, coll);
800
801 const int nm0 = m_base[0]->GetNumModes();
802 const int nm1 = m_base[1]->GetNumModes();
803 const int nm2 = m_base[2]->GetNumModes();
804
805 int mode0 = 0, mode1 = 0, mode2 = 0, cnt = 0;
806
807 bool found = false;
808 for (mode0 = 0; mode0 < nm0; ++mode0)
809 {
810 for (mode1 = 0; mode1 < nm1; ++mode1)
811 {
812 int maxpq = max(mode0, mode1);
813 for (mode2 = 0; mode2 < nm2 - maxpq; ++mode2, ++cnt)
814 {
815 if (cnt == mode)
816 {
817 found = true;
818 break;
819 }
820 }
821
822 if (found)
823 {
824 break;
825 }
826 }
827
828 if (found)
829 {
830 break;
831 }
832
833 for (int j = nm1; j < nm2; ++j)
834 {
835 int ijmax = max(mode0, j);
836 mode2 += nm2 - ijmax;
837 }
838 }
839
840 if (mode == 1 && m_base[0]->GetBasisType() == LibUtilities::eModified_A)
841 {
842 return StdExpansion::BaryEvaluateBasis<2>(coll[2], 1);
843 }
844 else
845 {
846 return StdExpansion::BaryEvaluateBasis<0>(coll[0], mode0) *
847 StdExpansion::BaryEvaluateBasis<1>(coll[1], mode1) *
848 StdExpansion::BaryEvaluateBasis<2>(coll[2], mode2);
849 }
850}

References Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::LocCoordToLocCollapsed(), and Nektar::StdRegions::StdExpansion::m_base.

◆ v_ReduceOrderCoeffs()

void Nektar::StdRegions::StdPyrExp::v_ReduceOrderCoeffs ( int  numMin,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2134 of file StdPyrExp.cpp.

2137{
2138 int nquad0 = m_base[0]->GetNumPoints();
2139 int nquad1 = m_base[1]->GetNumPoints();
2140 int nquad2 = m_base[2]->GetNumPoints();
2141 int nqtot = nquad0 * nquad1 * nquad2;
2142 int nmodes0 = m_base[0]->GetNumModes();
2143 int nmodes1 = m_base[1]->GetNumModes();
2144 int nmodes2 = m_base[2]->GetNumModes();
2145 int numMax = nmodes0;
2146
2147 Array<OneD, NekDouble> coeff(m_ncoeffs);
2148 Array<OneD, NekDouble> coeff_tmp1(m_ncoeffs, 0.0);
2149 Array<OneD, NekDouble> phys_tmp(nqtot, 0.0);
2150 Array<OneD, NekDouble> tmp, tmp2, tmp3, tmp4;
2151
2152 const LibUtilities::PointsKey Pkey0 = m_base[0]->GetPointsKey();
2153 const LibUtilities::PointsKey Pkey1 = m_base[1]->GetPointsKey();
2154 const LibUtilities::PointsKey Pkey2 = m_base[2]->GetPointsKey();
2155
2156 LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A, nmodes0, Pkey0);
2157 LibUtilities::BasisKey bortho1(LibUtilities::eOrtho_A, nmodes1, Pkey1);
2158 LibUtilities::BasisKey bortho2(LibUtilities::eOrthoPyr_C, nmodes2, Pkey2);
2159
2160 int cnt = 0;
2161 int u = 0;
2162 int i = 0;
2164
2166 bortho0, bortho1, bortho2);
2167
2168 BwdTrans(inarray, phys_tmp);
2169 OrthoPyrExp->FwdTrans(phys_tmp, coeff);
2170
2171 // filtering
2172 for (u = 0; u < numMin; ++u)
2173 {
2174 for (i = 0; i < numMin; ++i)
2175 {
2176
2177 int maxui = max(u, i);
2178 Vmath::Vcopy(numMin - maxui, tmp = coeff + cnt, 1,
2179 tmp2 = coeff_tmp1 + cnt, 1);
2180 cnt += nmodes2 - maxui;
2181 }
2182
2183 for (i = numMin; i < nmodes1; ++i)
2184 {
2185 int maxui = max(u, i);
2186 cnt += numMax - maxui;
2187 }
2188 }
2189
2190 OrthoPyrExp->BwdTrans(coeff_tmp1, phys_tmp);
2191 StdPyrExp::FwdTrans(phys_tmp, outarray);
2192}
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Backward transformation from coefficient space to physical space.
Definition: StdExpansion.h:424
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Forward transformation from physical space to coefficient space.
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:42
std::shared_ptr< StdPyrExp > StdPyrExpSharedPtr
Definition: StdPyrExp.h:214

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::StdRegions::StdExpansion::BwdTrans(), Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::eOrthoPyr_C, Nektar::StdRegions::StdExpansion::FwdTrans(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, and Vmath::Vcopy().

◆ v_StdPhysDeriv() [1/2]

void Nektar::StdRegions::StdPyrExp::v_StdPhysDeriv ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d0,
Array< OneD, NekDouble > &  out_d1,
Array< OneD, NekDouble > &  out_d2 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 188 of file StdPyrExp.cpp.

192{
193 StdPyrExp::v_PhysDeriv(inarray, out_d0, out_d1, out_d2);
194}

References v_PhysDeriv().

◆ v_StdPhysDeriv() [2/2]

void Nektar::StdRegions::StdPyrExp::v_StdPhysDeriv ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 196 of file StdPyrExp.cpp.

199{
200 StdPyrExp::v_PhysDeriv(dir, inarray, outarray);
201}

References v_PhysDeriv().

◆ v_SVVLaplacianFilter()

void Nektar::StdRegions::StdPyrExp::v_SVVLaplacianFilter ( Array< OneD, NekDouble > &  array,
const StdMatrixKey mkey 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PyrExp.

Definition at line 1995 of file StdPyrExp.cpp.

1997{
1998 // Generate an orthonogal expansion
1999 int qa = m_base[0]->GetNumPoints();
2000 int qb = m_base[1]->GetNumPoints();
2001 int qc = m_base[2]->GetNumPoints();
2002 int nmodes_a = m_base[0]->GetNumModes();
2003 int nmodes_b = m_base[1]->GetNumModes();
2004 int nmodes_c = m_base[2]->GetNumModes();
2005 // Declare orthogonal basis.
2006 LibUtilities::PointsKey pa(qa, m_base[0]->GetPointsType());
2007 LibUtilities::PointsKey pb(qb, m_base[1]->GetPointsType());
2008 LibUtilities::PointsKey pc(qc, m_base[2]->GetPointsType());
2009
2010 LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A, nmodes_a, pa);
2011 LibUtilities::BasisKey Bb(LibUtilities::eOrtho_A, nmodes_b, pb);
2012 LibUtilities::BasisKey Bc(LibUtilities::eOrthoPyr_C, nmodes_c, pc);
2013 StdPyrExp OrthoExp(Ba, Bb, Bc);
2014
2015 Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
2016 int i, j, k, cnt = 0;
2017
2018 // project onto modal space.
2019 OrthoExp.FwdTrans(array, orthocoeffs);
2020
2021 if (mkey.ConstFactorExists(eFactorSVVPowerKerDiffCoeff))
2022 {
2023 // Rodrigo's power kernel
2024 NekDouble cutoff = mkey.GetConstFactor(eFactorSVVCutoffRatio);
2025 NekDouble SvvDiffCoeff =
2026 mkey.GetConstFactor(eFactorSVVPowerKerDiffCoeff) *
2027 mkey.GetConstFactor(eFactorSVVDiffCoeff);
2028
2029 for (i = 0; i < nmodes_a; ++i)
2030 {
2031 for (j = 0; j < nmodes_b; ++j)
2032 {
2033 int maxij = max(i, j);
2034 NekDouble fac1 = std::max(
2035 pow((1.0 * i) / (nmodes_a - 1), cutoff * nmodes_a),
2036 pow((1.0 * j) / (nmodes_b - 1), cutoff * nmodes_b));
2037
2038 for (k = 0; k < nmodes_c - maxij; ++k)
2039 {
2040 NekDouble fac =
2041 std::max(fac1, pow((1.0 * k) / (nmodes_c - 1),
2042 cutoff * nmodes_c));
2043
2044 orthocoeffs[cnt] *= SvvDiffCoeff * fac;
2045 cnt++;
2046 }
2047 }
2048 }
2049 }
2050 else if (mkey.ConstFactorExists(
2051 eFactorSVVDGKerDiffCoeff)) // Rodrigo/Mansoor's DG Kernel
2052 {
2053 NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDGKerDiffCoeff) *
2054 mkey.GetConstFactor(eFactorSVVDiffCoeff);
2055
2056 int max_abc = max(nmodes_a - kSVVDGFiltermodesmin,
2057 nmodes_b - kSVVDGFiltermodesmin);
2058 max_abc = max(max_abc, nmodes_c - kSVVDGFiltermodesmin);
2059 // clamp max_abc
2060 max_abc = max(max_abc, 0);
2061 max_abc = min(max_abc, kSVVDGFiltermodesmax - kSVVDGFiltermodesmin);
2062
2063 for (i = 0; i < nmodes_a; ++i)
2064 {
2065 for (j = 0; j < nmodes_b; ++j)
2066 {
2067 int maxij = max(i, j);
2068
2069 for (k = 0; k < nmodes_c - maxij; ++k)
2070 {
2071 int maxijk = max(maxij, k);
2072 maxijk = min(maxijk, kSVVDGFiltermodesmax - 1);
2073
2074 orthocoeffs[cnt] *=
2075 SvvDiffCoeff * kSVVDGFilter[max_abc][maxijk];
2076 cnt++;
2077 }
2078 }
2079 }
2080 }
2081 else
2082 {
2083 // SVV filter paramaters (how much added diffusion relative
2084 // to physical one and fraction of modes from which you
2085 // start applying this added diffusion)
2086 //
2087 NekDouble SvvDiffCoeff =
2088 mkey.GetConstFactor(StdRegions::eFactorSVVDiffCoeff);
2089 NekDouble SVVCutOff =
2090 mkey.GetConstFactor(StdRegions::eFactorSVVCutoffRatio);
2091
2092 // Defining the cut of mode
2093 int cutoff_a = (int)(SVVCutOff * nmodes_a);
2094 int cutoff_b = (int)(SVVCutOff * nmodes_b);
2095 int cutoff_c = (int)(SVVCutOff * nmodes_c);
2096 // To avoid the fac[j] from blowing up
2097 NekDouble epsilon = 1;
2098
2099 int nmodes = min(min(nmodes_a, nmodes_b), nmodes_c);
2100 NekDouble cutoff = min(min(cutoff_a, cutoff_b), cutoff_c);
2101
2102 for (i = 0; i < nmodes_a; ++i) // P
2103 {
2104 for (j = 0; j < nmodes_b; ++j) // Q
2105 {
2106 int maxij = max(i, j);
2107 for (k = 0; k < nmodes_c - maxij; ++k) // R
2108 {
2109 if (j + k >= cutoff || i + k >= cutoff)
2110 {
2111 orthocoeffs[cnt] *=
2112 (SvvDiffCoeff *
2113 exp(-(i + k - nmodes) * (i + k - nmodes) /
2114 ((NekDouble)((i + k - cutoff + epsilon) *
2115 (i + k - cutoff + epsilon)))) *
2116 exp(-(j - nmodes) * (j - nmodes) /
2117 ((NekDouble)((j - cutoff + epsilon) *
2118 (j - cutoff + epsilon)))));
2119 }
2120 else
2121 {
2122 orthocoeffs[cnt] *= 0.0;
2123 }
2124 cnt++;
2125 }
2126 }
2127 }
2128 }
2129
2130 // backward transform to physical space
2131 OrthoExp.BwdTrans(orthocoeffs, array);
2132}
const int kSVVDGFiltermodesmin
Definition: StdRegions.hpp:500
const int kSVVDGFiltermodesmax
Definition: StdRegions.hpp:501
const NekDouble kSVVDGFilter[9][11]
Definition: StdRegions.hpp:503

References Nektar::StdRegions::StdExpansion::BwdTrans(), Nektar::StdRegions::StdMatrixKey::ConstFactorExists(), Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::eFactorSVVDGKerDiffCoeff, Nektar::StdRegions::eFactorSVVDiffCoeff, Nektar::StdRegions::eFactorSVVPowerKerDiffCoeff, Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::eOrthoPyr_C, Nektar::StdRegions::StdExpansion::FwdTrans(), Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::StdRegions::kSVVDGFilter, Nektar::StdRegions::kSVVDGFiltermodesmax, Nektar::StdRegions::kSVVDGFiltermodesmin, and Nektar::StdRegions::StdExpansion::m_base.