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

Class representing a prismatic element in reference space. More...

#include <StdPrismExp.h>

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

Public Member Functions

 StdPrismExp ()
 
 StdPrismExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
 
 StdPrismExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc, NekDouble *coeffs, NekDouble *phys)
 
 StdPrismExp (const StdPrismExp &T)
 
 ~StdPrismExp ()
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion3D
 StdExpansion3D ()
 
 StdExpansion3D (int numcoeffs, const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
 
 StdExpansion3D (const StdExpansion3D &T)
 
virtual ~StdExpansion3D ()
 
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 (void) const
 
std::shared_ptr< StdExpansionGetLinStdExp (void) const
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion ()
 
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 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, 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 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

virtual void v_PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
 Calculate the derivative of the physical points. More...
 
virtual void v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Calculate the derivative of the physical points in a given direction. More...
 
virtual void v_StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
 
virtual void v_StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
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)
 
virtual void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Forward transform from physical quadrature space stored in inarray and evaluate the expansion coefficients and store in outarray. More...
 
virtual void v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Calculate the inner product of inarray with respect to the basis B=base0*base1*base2 and put into outarray: More...
 
virtual void v_IProductWRTBase_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
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)
 
virtual void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Inner product of inarray over region with respect to the object's default expansion basis; output in outarray. More...
 
virtual void v_IProductWRTDerivBase_MatOp (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 
virtual void v_LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi)
 
virtual void v_GetCoords (Array< OneD, NekDouble > &xi_x, Array< OneD, NekDouble > &xi_y, Array< OneD, NekDouble > &xi_z)
 
virtual void v_FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 
NekDouble v_PhysEvaluateBasis (const Array< OneD, const NekDouble > &coords, int mode) final
 
virtual void v_GetTraceNumModes (const int fid, int &numModes0, int &numModes1, Orientation faceOrient=eDir1FwdDir1_Dir2FwdDir2)
 
virtual int v_GetNverts () const
 
virtual int v_GetNedges () const
 
virtual int v_GetNtraces () const
 
virtual LibUtilities::ShapeType v_DetShapeType () const
 Return Shape of region, using ShapeType enum list; i.e. prism. More...
 
virtual int v_NumBndryCoeffs () const
 
virtual int v_NumDGBndryCoeffs () const
 
virtual int v_GetTraceNcoeffs (const int i) const
 
virtual int v_GetTotalTraceIntNcoeffs () const
 
virtual int v_GetTraceIntNcoeffs (const int i) const
 
virtual int v_GetTraceNumPoints (const int i) const
 
virtual int v_GetEdgeNcoeffs (const int i) const
 
virtual const LibUtilities::BasisKey v_GetTraceBasisKey (const int i, const int k) const
 
virtual LibUtilities::PointsKey v_GetTracePointsKey (const int i, const int j) const
 
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
virtual bool v_IsBoundaryInteriorExpansion ()
 
virtual int v_GetVertexMap (int localVertexId, bool useCoeffPacking=false)
 
virtual void v_GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
virtual void v_GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
virtual void v_GetTraceCoeffMap (const unsigned int fid, Array< OneD, unsigned int > &maparray)
 
virtual void v_GetElmtTraceToTraceMap (const unsigned int fid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation faceOrient, int P, int Q)
 
virtual void v_GetEdgeInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2)
 
virtual void v_GetTraceInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2)
 
virtual DNekMatSharedPtr v_GenMatrix (const StdMatrixKey &mkey)
 
virtual DNekMatSharedPtr v_CreateStdMatrix (const StdMatrixKey &mkey)
 
virtual void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void v_SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
virtual void v_ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion3D
virtual NekDouble v_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...
 
virtual NekDouble v_PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 
virtual void v_LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
 
virtual void v_HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
 
virtual NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray)
 Integrates the specified function over the domain. More...
 
virtual void v_GetTraceToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient, int P, int Q)
 
virtual void v_GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 
- 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_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
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 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)
 
template<int DIR>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals)
 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)
 

Private Member Functions

int GetMode (int I, int J, int K)
 Compute the local 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

Class representing a prismatic element in reference space.

Definition at line 48 of file StdPrismExp.h.

Constructor & Destructor Documentation

◆ StdPrismExp() [1/4]

Nektar::StdRegions::StdPrismExp::StdPrismExp ( )

Definition at line 47 of file StdPrismExp.cpp.

49 {
50 }

◆ StdPrismExp() [2/4]

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

Definition at line 52 of file StdPrismExp.cpp.

56  Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
57  3, Ba, Bb, Bc),
59  Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
60  Ba, Bb, Bc)
61 {
62  ASSERTL0(Ba.GetNumModes() <= Bc.GetNumModes(),
63  "order in 'a' direction is higher than order in 'c' direction");
64 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
StdExpansion()
Default Constructor.
int getNumberOfCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:282

References ASSERTL0, and Nektar::LibUtilities::BasisKey::GetNumModes().

◆ StdPrismExp() [3/4]

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

◆ StdPrismExp() [4/4]

Nektar::StdRegions::StdPrismExp::StdPrismExp ( const StdPrismExp T)

Definition at line 66 of file StdPrismExp.cpp.

68 {
69 }

◆ ~StdPrismExp()

Nektar::StdRegions::StdPrismExp::~StdPrismExp ( )

Definition at line 72 of file StdPrismExp.cpp.

73 {
74 }

Member Function Documentation

◆ GetMode()

int Nektar::StdRegions::StdPrismExp::GetMode ( int  p,
int  q,
int  r 
)
private

Compute the local 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). For example, with P=1, Q=2, R=3, the indexing inside each q-r plane (with r increasing upwards and q to the right) is:

p = 0: p = 1:

3 7 11 2 6 10 14 17 20 1 5 9 13 16 19 0 4 8 12 15 18

Note that in this element, we must have that \( P <= R \).

Definition at line 1958 of file StdPrismExp.cpp.

1959 {
1960  int Q = m_base[1]->GetNumModes() - 1;
1961  int R = m_base[2]->GetNumModes() - 1;
1962 
1963  return r + // Skip along stacks (r-direction)
1964  q * (R + 1 - p) + // Skip along columns (q-direction)
1965  (Q + 1) * (p * R + 1 -
1966  (p - 2) * (p - 1) / 2); // Skip along rows (p-direction)
1967 }
Array< OneD, LibUtilities::BasisSharedPtr > m_base

References Nektar::StdRegions::StdExpansion::m_base, and CellMLToNektar.cellml_metadata::p.

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

◆ v_BwdTrans()

void Nektar::StdRegions::StdPrismExp::v_BwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual
Note
'r' (base[2]) runs fastest in this element.

Perform backwards transformation 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})\)

In the prism this expansion becomes:

\( 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_{pr}^b (\xi_{3k}) \rbrace} \rbrace}. \)

And sumfactorizing step of the form is as:\

\( f_{pr} (\xi_{3k}) = \sum_{r=0}^{Q_z} \hat u_{pqr} \psi_{pr}^b (\xi_{3k}),\\ g_{p} (\xi_{2j}, \xi_{3k}) = \sum_{r=0}^{Q_y} \psi_{p}^a (\xi_{2j}) f_{pr} (\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.

Reimplemented in Nektar::StdRegions::StdNodalPrismExp.

Definition at line 239 of file StdPrismExp.cpp.

241 {
244  "Basis[1] is not a general tensor type");
245 
248  "Basis[2] is not a general tensor type");
249 
250  if (m_base[0]->Collocation() && m_base[1]->Collocation() &&
251  m_base[2]->Collocation())
252  {
254  m_base[2]->GetNumPoints(),
255  inarray, 1, outarray, 1);
256  }
257  else
258  {
259  StdPrismExp::v_BwdTrans_SumFac(inarray, outarray);
260  }
261 }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:163
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:226
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:51
@ eModified_C
Principle Modified Functions .
Definition: BasisType.h:52
@ eOrtho_C
Principle Orthogonal Functions .
Definition: BasisType.h:48
@ eOrtho_B
Principle Orthogonal Functions .
Definition: BasisType.h:46
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255

References ASSERTL1, Nektar::LibUtilities::eModified_B, Nektar::LibUtilities::eModified_C, Nektar::LibUtilities::eOrtho_B, Nektar::LibUtilities::eOrtho_C, Nektar::StdRegions::StdExpansion::GetBasisType(), 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::StdPrismExp::v_BwdTrans_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalPrismExp.

Definition at line 263 of file StdPrismExp.cpp.

265 {
266  int nquad1 = m_base[1]->GetNumPoints();
267  int nquad2 = m_base[2]->GetNumPoints();
268  int order0 = m_base[0]->GetNumModes();
269  int order1 = m_base[1]->GetNumModes();
270 
271  Array<OneD, NekDouble> wsp(nquad2 * order1 * order0 +
272  nquad1 * nquad2 * order0);
273 
274  BwdTrans_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetBdata(),
275  m_base[2]->GetBdata(), inarray, outarray, wsp, true,
276  true, true);
277 }
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)

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

Referenced by v_BwdTrans(), and Nektar::StdRegions::StdNodalPrismExp::v_BwdTrans_SumFac().

◆ v_BwdTrans_SumFacKernel()

void Nektar::StdRegions::StdPrismExp::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 
)
protectedvirtual

Implements Nektar::StdRegions::StdExpansion3D.

Definition at line 279 of file StdPrismExp.cpp.

286 {
287  boost::ignore_unused(doCheckCollDir0, doCheckCollDir1, doCheckCollDir2);
288 
289  int i, mode;
290  int nquad0 = m_base[0]->GetNumPoints();
291  int nquad1 = m_base[1]->GetNumPoints();
292  int nquad2 = m_base[2]->GetNumPoints();
293  int nummodes0 = m_base[0]->GetNumModes();
294  int nummodes1 = m_base[1]->GetNumModes();
295  int nummodes2 = m_base[2]->GetNumModes();
296  Array<OneD, NekDouble> tmp0 = wsp;
297  Array<OneD, NekDouble> tmp1 = tmp0 + nquad2 * nummodes1 * nummodes0;
298 
299  for (i = mode = 0; i < nummodes0; ++i)
300  {
301  Blas::Dgemm('N', 'N', nquad2, nummodes1, nummodes2 - i, 1.0,
302  base2.get() + mode * nquad2, nquad2,
303  inarray.get() + mode * nummodes1, nummodes2 - i, 0.0,
304  tmp0.get() + i * nquad2 * nummodes1, nquad2);
305  mode += nummodes2 - i;
306  }
307 
309  {
310  for (i = 0; i < nummodes1; i++)
311  {
312  Blas::Daxpy(nquad2, inarray[1 + i * nummodes2],
313  base2.get() + nquad2, 1,
314  tmp0.get() + nquad2 * (nummodes1 + i), 1);
315  }
316  }
317 
318  for (i = 0; i < nummodes0; i++)
319  {
320  Blas::Dgemm('N', 'T', nquad1, nquad2, nummodes1, 1.0, base1.get(),
321  nquad1, tmp0.get() + i * nquad2 * nummodes1, nquad2, 0.0,
322  tmp1.get() + i * nquad2 * nquad1, nquad1);
323  }
324 
325  Blas::Dgemm('N', 'T', nquad0, nquad2 * nquad1, nummodes0, 1.0, base0.get(),
326  nquad0, tmp1.get(), nquad2 * nquad1, 0.0, outarray.get(),
327  nquad0);
328 }
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:368
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:154
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:50

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

◆ v_CalcNumberOfCoefficients()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 983 of file StdPrismExp.cpp.

985 {
987  nummodes[modes_offset], nummodes[modes_offset + 1],
988  nummodes[modes_offset + 2]);
989 
990  modes_offset += 3;
991  return nmodes;
992 }

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

◆ v_CreateStdMatrix()

DNekMatSharedPtr Nektar::StdRegions::StdPrismExp::v_CreateStdMatrix ( const StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PrismExp, and Nektar::StdRegions::StdNodalPrismExp.

Definition at line 1935 of file StdPrismExp.cpp.

1936 {
1937  return v_GenMatrix(mkey);
1938 }
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey)

References v_GenMatrix().

◆ v_DetShapeType()

LibUtilities::ShapeType Nektar::StdRegions::StdPrismExp::v_DetShapeType ( ) const
protectedvirtual

Return Shape of region, using ShapeType enum list; i.e. prism.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 816 of file StdPrismExp.cpp.

817 {
818  return LibUtilities::ePrism;
819 }

References Nektar::LibUtilities::ePrism.

◆ v_FillMode()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalPrismExp.

Definition at line 699 of file StdPrismExp.cpp.

700 {
701  Array<OneD, NekDouble> tmp(m_ncoeffs, 0.0);
702  tmp[mode] = 1.0;
703  StdPrismExp::v_BwdTrans(tmp, outarray);
704 }
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)

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

Referenced by Nektar::StdRegions::StdNodalPrismExp::GenNBasisTransMatrix().

◆ v_FwdTrans()

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

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::StdRegions::StdNodalPrismExp, and Nektar::LocalRegions::PrismExp.

Definition at line 341 of file StdPrismExp.cpp.

343 {
344  v_IProductWRTBase(inarray, outarray);
345 
346  // Get Mass matrix inverse
347  StdMatrixKey masskey(eInvMass, DetShapeType(), *this);
348  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
349 
350  // copy inarray in case inarray == outarray
351  DNekVec in(m_ncoeffs, outarray);
352  DNekVec out(m_ncoeffs, outarray, eWrapper);
353 
354  out = (*matsys) * in;
355 }
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:611
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:375
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Calculate the inner product of inarray with respect to the basis B=base0*base1*base2 and put into out...
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::StdPrismExp::v_GenMatrix ( const StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PrismExp, and Nektar::StdRegions::StdNodalPrismExp.

Definition at line 1848 of file StdPrismExp.cpp.

1849 {
1850 
1851  MatrixType mtype = mkey.GetMatrixType();
1852 
1853  DNekMatSharedPtr Mat;
1854 
1855  switch (mtype)
1856  {
1858  {
1859  int nq0 = m_base[0]->GetNumPoints();
1860  int nq1 = m_base[1]->GetNumPoints();
1861  int nq2 = m_base[2]->GetNumPoints();
1862  int nq;
1863 
1864  // take definition from key
1865  if (mkey.ConstFactorExists(eFactorConst))
1866  {
1867  nq = (int)mkey.GetConstFactor(eFactorConst);
1868  }
1869  else
1870  {
1871  nq = max(nq0, max(nq1, nq2));
1872  }
1873 
1874  int neq =
1876  Array<OneD, Array<OneD, NekDouble>> coords(neq);
1877  Array<OneD, NekDouble> coll(3);
1878  Array<OneD, DNekMatSharedPtr> I(3);
1879  Array<OneD, NekDouble> tmp(nq0);
1880 
1881  Mat =
1882  MemoryManager<DNekMat>::AllocateSharedPtr(neq, nq0 * nq1 * nq2);
1883  int cnt = 0;
1884  for (int i = 0; i < nq; ++i)
1885  {
1886  for (int j = 0; j < nq; ++j)
1887  {
1888  for (int k = 0; k < nq - i; ++k, ++cnt)
1889  {
1890  coords[cnt] = Array<OneD, NekDouble>(3);
1891  coords[cnt][0] = -1.0 + 2 * k / (NekDouble)(nq - 1);
1892  coords[cnt][1] = -1.0 + 2 * j / (NekDouble)(nq - 1);
1893  coords[cnt][2] = -1.0 + 2 * i / (NekDouble)(nq - 1);
1894  }
1895  }
1896  }
1897 
1898  for (int i = 0; i < neq; ++i)
1899  {
1900  LocCoordToLocCollapsed(coords[i], coll);
1901 
1902  I[0] = m_base[0]->GetI(coll);
1903  I[1] = m_base[1]->GetI(coll + 1);
1904  I[2] = m_base[2]->GetI(coll + 2);
1905 
1906  // interpolate first coordinate direction
1907  NekDouble fac;
1908  for (int k = 0; k < nq2; ++k)
1909  {
1910  for (int j = 0; j < nq1; ++j)
1911  {
1912 
1913  fac = (I[1]->GetPtr())[j] * (I[2]->GetPtr())[k];
1914  Vmath::Smul(nq0, fac, I[0]->GetPtr(), 1, tmp, 1);
1915 
1916  Vmath::Vcopy(nq0, &tmp[0], 1,
1917  Mat->GetRawPtr() + k * nq0 * nq1 * neq +
1918  j * nq0 * neq + i,
1919  neq);
1920  }
1921  }
1922  }
1923  }
1924  break;
1925  default:
1926  {
1928  }
1929  break;
1930  }
1931 
1932  return Mat;
1933 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
Definition: StdExpansion.h:974
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
double NekDouble
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.cpp:248

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::StdRegions::StdMatrixKey::ConstFactorExists(), Nektar::StdRegions::StdExpansion::CreateGeneralMatrix(), Nektar::StdRegions::eFactorConst, Nektar::StdRegions::ePhysInterpToEquiSpaced, Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdMatrixKey::GetMatrixType(), Nektar::LibUtilities::StdPrismData::getNumberOfCoefficients(), Nektar::StdRegions::StdExpansion::LocCoordToLocCollapsed(), Nektar::StdRegions::StdExpansion::m_base, Vmath::Smul(), and Vmath::Vcopy().

Referenced by v_CreateStdMatrix().

◆ v_GetBoundaryMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalPrismExp.

Definition at line 1108 of file StdPrismExp.cpp.

1109 {
1112  "BasisType is not a boundary interior form");
1115  "BasisType is not a boundary interior form");
1118  "BasisType is not a boundary interior form");
1119 
1120  int P = m_base[0]->GetNumModes() - 1, p;
1121  int Q = m_base[1]->GetNumModes() - 1, q;
1122  int R = m_base[2]->GetNumModes() - 1, r;
1123  int idx = 0;
1124 
1125  int nBnd = NumBndryCoeffs();
1126 
1127  if (maparray.size() != nBnd)
1128  {
1129  maparray = Array<OneD, unsigned int>(nBnd);
1130  }
1131 
1132  // Loop over all boundary modes (in ascending order).
1133  for (p = 0; p <= P; ++p)
1134  {
1135  // First two q-r planes are entirely boundary modes.
1136  if (p <= 1)
1137  {
1138  for (q = 0; q <= Q; ++q)
1139  {
1140  for (r = 0; r <= R - p; ++r)
1141  {
1142  maparray[idx++] = GetMode(p, q, r);
1143  }
1144  }
1145  }
1146  else
1147  {
1148  // Remaining q-r planes contain boundary modes on the two
1149  // left-hand sides and bottom edge.
1150  for (q = 0; q <= Q; ++q)
1151  {
1152  if (q <= 1)
1153  {
1154  for (r = 0; r <= R - p; ++r)
1155  {
1156  maparray[idx++] = GetMode(p, q, r);
1157  }
1158  }
1159  else
1160  {
1161  maparray[idx++] = GetMode(p, q, 0);
1162  }
1163  }
1164  }
1165  }
1166 }
int GetMode(int I, int J, int K)
Compute the local mode number in the expansion for a particular tensorial combination.
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:58

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

◆ v_GetCoords()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PrismExp.

Definition at line 672 of file StdPrismExp.cpp.

675 {
676  Array<OneD, const NekDouble> etaBar_x = m_base[0]->GetZ();
677  Array<OneD, const NekDouble> eta_y = m_base[1]->GetZ();
678  Array<OneD, const NekDouble> eta_z = m_base[2]->GetZ();
679  int Qx = GetNumPoints(0);
680  int Qy = GetNumPoints(1);
681  int Qz = GetNumPoints(2);
682 
683  // Convert collapsed coordinates into cartesian coordinates: eta --> xi
684  for (int k = 0; k < Qz; ++k)
685  {
686  for (int j = 0; j < Qy; ++j)
687  {
688  for (int i = 0; i < Qx; ++i)
689  {
690  int s = i + Qx * (j + Qy * k);
691  xi_x[s] = (1.0 - eta_z[k]) * (1.0 + etaBar_x[i]) / 2.0 - 1.0;
692  xi_y[s] = eta_y[j];
693  xi_z[s] = eta_z[k];
694  }
695  }
696  }
697 }

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

◆ v_GetEdgeInteriorToElementMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 1545 of file StdPrismExp.cpp.

1548 {
1549  int i;
1550  bool signChange;
1551  const int P = m_base[0]->GetNumModes() - 1;
1552  const int Q = m_base[1]->GetNumModes() - 1;
1553  const int R = m_base[2]->GetNumModes() - 1;
1554  const int nEdgeIntCoeffs = v_GetEdgeNcoeffs(eid) - 2;
1555 
1556  if (maparray.size() != nEdgeIntCoeffs)
1557  {
1558  maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
1559  }
1560 
1561  if (signarray.size() != nEdgeIntCoeffs)
1562  {
1563  signarray = Array<OneD, int>(nEdgeIntCoeffs, 1);
1564  }
1565  else
1566  {
1567  fill(signarray.get(), signarray.get() + nEdgeIntCoeffs, 1);
1568  }
1569 
1570  // If edge is oriented backwards, change sign of modes which have
1571  // degree 2n+1, n >= 1.
1572  signChange = edgeOrient == eBackwards;
1573 
1574  switch (eid)
1575  {
1576  case 0:
1577  for (i = 2; i <= P; ++i)
1578  {
1579  maparray[i - 2] = GetMode(i, 0, 0);
1580  }
1581  break;
1582 
1583  case 1:
1584  for (i = 2; i <= Q; ++i)
1585  {
1586  maparray[i - 2] = GetMode(1, i, 0);
1587  }
1588  break;
1589 
1590  case 2:
1591  // Base quad; reverse direction.
1592  // signChange = !signChange;
1593  for (i = 2; i <= P; ++i)
1594  {
1595  maparray[i - 2] = GetMode(i, 1, 0);
1596  }
1597  break;
1598 
1599  case 3:
1600  // Base quad; reverse direction.
1601  // signChange = !signChange;
1602  for (i = 2; i <= Q; ++i)
1603  {
1604  maparray[i - 2] = GetMode(0, i, 0);
1605  }
1606  break;
1607 
1608  case 4:
1609  for (i = 2; i <= R; ++i)
1610  {
1611  maparray[i - 2] = GetMode(0, 0, i);
1612  }
1613  break;
1614 
1615  case 5:
1616  for (i = 1; i <= R - 1; ++i)
1617  {
1618  maparray[i - 1] = GetMode(1, 0, i);
1619  }
1620  break;
1621 
1622  case 6:
1623  for (i = 1; i <= R - 1; ++i)
1624  {
1625  maparray[i - 1] = GetMode(1, 1, i);
1626  }
1627  break;
1628 
1629  case 7:
1630  for (i = 2; i <= R; ++i)
1631  {
1632  maparray[i - 2] = GetMode(0, 1, i);
1633  }
1634  break;
1635 
1636  case 8:
1637  for (i = 2; i <= Q; ++i)
1638  {
1639  maparray[i - 2] = GetMode(0, i, 1);
1640  }
1641  break;
1642 
1643  default:
1644  ASSERTL0(false, "Edge not defined.");
1645  break;
1646  }
1647 
1648  if (signChange)
1649  {
1650  for (i = 1; i < nEdgeIntCoeffs; i += 2)
1651  {
1652  signarray[i] = -1;
1653  }
1654  }
1655 }
virtual int v_GetEdgeNcoeffs(const int i) const

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

◆ v_GetEdgeNcoeffs()

int Nektar::StdRegions::StdPrismExp::v_GetEdgeNcoeffs ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 775 of file StdPrismExp.cpp.

776 {
777  ASSERTL2(i >= 0 && i <= 8, "edge id is out of range");
778 
779  if (i == 0 || i == 2)
780  {
781  return GetBasisNumModes(0);
782  }
783  else if (i == 1 || i == 3 || i == 8)
784  {
785  return GetBasisNumModes(1);
786  }
787  else
788  {
789  return GetBasisNumModes(2);
790  }
791 }
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:272
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:176

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

Referenced by v_GetEdgeInteriorToElementMap().

◆ v_GetElmtTraceToTraceMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1269 of file StdPrismExp.cpp.

1273 {
1275  "Method only implemented if BasisType is identical"
1276  "in x and y directions");
1279  "Method only implemented for Modified_A BasisType"
1280  "(x and y direction) and Modified_B BasisType (z "
1281  "direction)");
1282 
1283  int i, j, k, p, r, nFaceCoeffs, idx = 0;
1284  int nummodesA = 0, nummodesB = 0;
1285 
1286  switch (fid)
1287  {
1288  case 0:
1289  nummodesA = m_base[0]->GetNumModes();
1290  nummodesB = m_base[1]->GetNumModes();
1291  break;
1292  case 1:
1293  case 3:
1294  nummodesA = m_base[0]->GetNumModes();
1295  nummodesB = m_base[2]->GetNumModes();
1296  break;
1297  case 2:
1298  case 4:
1299  nummodesA = m_base[1]->GetNumModes();
1300  nummodesB = m_base[2]->GetNumModes();
1301  break;
1302  default:
1303  ASSERTL0(false, "fid must be between 0 and 4");
1304  }
1305 
1306  if (P == -1)
1307  {
1308  P = nummodesA;
1309  Q = nummodesB;
1310  nFaceCoeffs = GetTraceNcoeffs(fid);
1311  }
1312  else if (fid == 1 || fid == 3)
1313  {
1314  nFaceCoeffs = P * (2 * Q - P + 1) / 2;
1315  }
1316  else
1317  {
1318  nFaceCoeffs = P * Q;
1319  }
1320 
1321  // Allocate the map array and sign array; set sign array to ones (+)
1322  if (maparray.size() != nFaceCoeffs)
1323  {
1324  maparray = Array<OneD, unsigned int>(nFaceCoeffs);
1325  }
1326 
1327  if (signarray.size() != nFaceCoeffs)
1328  {
1329  signarray = Array<OneD, int>(nFaceCoeffs, 1);
1330  }
1331  else
1332  {
1333  fill(signarray.get(), signarray.get() + nFaceCoeffs, 1);
1334  }
1335 
1336  int minPA = min(nummodesA, P);
1337  int minQB = min(nummodesB, Q);
1338  // triangular faces
1339  if (fid == 1 || fid == 3)
1340  {
1341  // zero signmap and set maparray to zero if elemental
1342  // modes are not as large as face modesl
1343  idx = 0;
1344  int cnt = 0;
1345 
1346  for (j = 0; j < minPA; ++j)
1347  {
1348  // set maparray
1349  for (k = 0; k < minQB - j; ++k, ++cnt)
1350  {
1351  maparray[idx++] = cnt;
1352  }
1353 
1354  cnt += nummodesB - minQB;
1355 
1356  // idx += nummodesB-j;
1357  for (k = nummodesB - j; k < Q - j; ++k)
1358  {
1359  signarray[idx] = 0.0;
1360  maparray[idx++] = maparray[0];
1361  }
1362  }
1363 #if 0 // no required?
1364  for (j = minPA; j < nummodesA; ++j)
1365  {
1366  // set maparray
1367  for (k = 0; k < minQB-j; ++k, ++cnt)
1368  {
1369  maparray[idx++] = cnt;
1370  }
1371 
1372  cnt += nummodesB-minQB;
1373 
1374  //idx += nummodesB-j;
1375  for (k = nummodesB-j; k < Q-j; ++k)
1376  {
1377  signarray[idx] = 0.0;
1378  maparray[idx++] = maparray[0];
1379  }
1380  }
1381 #endif
1382  for (j = nummodesA; j < P; ++j)
1383  {
1384  for (k = 0; k < Q - j; ++k)
1385  {
1386  signarray[idx] = 0.0;
1387  maparray[idx++] = maparray[0];
1388  }
1389  }
1390 
1391  // Triangles only have one possible orientation (base
1392  // direction reversed); swap edge modes.
1393  if (faceOrient == eDir1BwdDir1_Dir2FwdDir2)
1394  {
1395  idx = 0;
1396  for (p = 0; p < P; ++p)
1397  {
1398  for (r = 0; r < Q - p; ++r, idx++)
1399  {
1400  if (p > 1)
1401  {
1402  signarray[idx] = p % 2 ? -1 : 1;
1403  }
1404  }
1405  }
1406 
1407  swap(maparray[0], maparray[Q]);
1408  for (i = 1; i < Q - 1; ++i)
1409  {
1410  swap(maparray[i + 1], maparray[Q + i]);
1411  }
1412  }
1413  }
1414  else
1415  {
1416  // Set up an array indexing for quads, since the
1417  // ordering may need to be transposed.
1418  Array<OneD, int> arrayindx(nFaceCoeffs, -1);
1419 
1420  for (i = 0; i < Q; i++)
1421  {
1422  for (j = 0; j < P; j++)
1423  {
1424  if (faceOrient < eDir1FwdDir2_Dir2FwdDir1)
1425  {
1426  arrayindx[i * P + j] = i * P + j;
1427  }
1428  else
1429  {
1430  arrayindx[i * P + j] = j * Q + i;
1431  }
1432  }
1433  }
1434 
1435  // zero signmap and set maparray to zero if elemental
1436  // modes are not as large as face modesl
1437  for (j = 0; j < P; ++j)
1438  {
1439  // set up default maparray
1440  for (k = 0; k < Q; k++)
1441  {
1442  maparray[arrayindx[j + k * P]] = j + k * nummodesA;
1443  }
1444 
1445  for (k = nummodesB; k < Q; ++k)
1446  {
1447  signarray[arrayindx[j + k * P]] = 0.0;
1448  maparray[arrayindx[j + k * P]] = maparray[0];
1449  }
1450  }
1451 
1452  for (j = nummodesA; j < P; ++j)
1453  {
1454  for (k = 0; k < Q; ++k)
1455  {
1456  signarray[arrayindx[j + k * P]] = 0.0;
1457  maparray[arrayindx[j + k * P]] = maparray[0];
1458  }
1459  }
1460 
1461  // The code below is exactly the same as that taken from
1462  // StdHexExp and reverses the 'b' and 'a' directions as
1463  // appropriate (1st and 2nd if statements respectively) in
1464  // quadrilateral faces.
1465  if (faceOrient == eDir1FwdDir1_Dir2BwdDir2 ||
1466  faceOrient == eDir1BwdDir1_Dir2BwdDir2 ||
1467  faceOrient == eDir1BwdDir2_Dir2FwdDir1 ||
1468  faceOrient == eDir1BwdDir2_Dir2BwdDir1)
1469  {
1470  if (faceOrient < eDir1FwdDir2_Dir2FwdDir1)
1471  {
1472  for (i = 3; i < Q; i += 2)
1473  {
1474  for (j = 0; j < P; j++)
1475  {
1476  signarray[arrayindx[i * P + j]] *= -1;
1477  }
1478  }
1479 
1480  for (i = 0; i < P; i++)
1481  {
1482  swap(maparray[i], maparray[i + P]);
1483  swap(signarray[i], signarray[i + P]);
1484  }
1485  }
1486  else
1487  {
1488  for (i = 0; i < Q; i++)
1489  {
1490  for (j = 3; j < P; j += 2)
1491  {
1492  signarray[arrayindx[i * P + j]] *= -1;
1493  }
1494  }
1495 
1496  for (i = 0; i < Q; i++)
1497  {
1498  swap(maparray[i], maparray[i + Q]);
1499  swap(signarray[i], signarray[i + Q]);
1500  }
1501  }
1502  }
1503 
1504  if (faceOrient == eDir1BwdDir1_Dir2FwdDir2 ||
1505  faceOrient == eDir1BwdDir1_Dir2BwdDir2 ||
1506  faceOrient == eDir1FwdDir2_Dir2BwdDir1 ||
1507  faceOrient == eDir1BwdDir2_Dir2BwdDir1)
1508  {
1509  if (faceOrient < eDir1FwdDir2_Dir2FwdDir1)
1510  {
1511  for (i = 0; i < Q; i++)
1512  {
1513  for (j = 3; j < P; j += 2)
1514  {
1515  signarray[arrayindx[i * P + j]] *= -1;
1516  }
1517  }
1518 
1519  for (i = 0; i < Q; i++)
1520  {
1521  swap(maparray[i * P], maparray[i * P + 1]);
1522  swap(signarray[i * P], signarray[i * P + 1]);
1523  }
1524  }
1525  else
1526  {
1527  for (i = 3; i < Q; i += 2)
1528  {
1529  for (j = 0; j < P; j++)
1530  {
1531  signarray[arrayindx[i * P + j]] *= -1;
1532  }
1533  }
1534 
1535  for (i = 0; i < P; i++)
1536  {
1537  swap(maparray[i * Q], maparray[i * Q + 1]);
1538  swap(signarray[i * Q], signarray[i * Q + 1]);
1539  }
1540  }
1541  }
1542  }
1543 }
int GetTraceNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th trace.
Definition: StdExpansion.h:269

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::eModified_B, 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::StdPrismExp::v_GetInteriorMap ( Array< OneD, unsigned int > &  outarray)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalPrismExp.

Definition at line 1070 of file StdPrismExp.cpp.

1071 {
1074  "BasisType is not a boundary interior form");
1077  "BasisType is not a boundary interior form");
1080  "BasisType is not a boundary interior form");
1081 
1082  int P = m_base[0]->GetNumModes() - 1, p;
1083  int Q = m_base[1]->GetNumModes() - 1, q;
1084  int R = m_base[2]->GetNumModes() - 1, r;
1085 
1086  int nIntCoeffs = m_ncoeffs - NumBndryCoeffs();
1087 
1088  if (outarray.size() != nIntCoeffs)
1089  {
1090  outarray = Array<OneD, unsigned int>(nIntCoeffs);
1091  }
1092 
1093  int idx = 0;
1094 
1095  // Loop over all interior modes.
1096  for (p = 2; p <= P; ++p)
1097  {
1098  for (q = 2; q <= Q; ++q)
1099  {
1100  for (r = 1; r <= R - p; ++r)
1101  {
1102  outarray[idx++] = GetMode(p, q, r);
1103  }
1104  }
1105  }
1106 }

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

◆ v_GetNedges()

int Nektar::StdRegions::StdPrismExp::v_GetNedges ( void  ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 802 of file StdPrismExp.cpp.

803 {
804  return 9;
805 }

◆ v_GetNtraces()

int Nektar::StdRegions::StdPrismExp::v_GetNtraces ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 807 of file StdPrismExp.cpp.

808 {
809  return 5;
810 }

◆ v_GetNverts()

int Nektar::StdRegions::StdPrismExp::v_GetNverts ( ) const
protectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 797 of file StdPrismExp.cpp.

798 {
799  return 6;
800 }

◆ v_GetTotalTraceIntNcoeffs()

int Nektar::StdRegions::StdPrismExp::v_GetTotalTraceIntNcoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 901 of file StdPrismExp.cpp.

902 {
903  int Pi = GetBasisNumModes(0) - 2;
904  int Qi = GetBasisNumModes(1) - 2;
905  int Ri = GetBasisNumModes(2) - 2;
906 
907  return Pi * Qi + Pi * (2 * Ri - Pi - 1) + 2 * Qi * Ri;
908 }

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

◆ v_GetTraceBasisKey()

const LibUtilities::BasisKey Nektar::StdRegions::StdPrismExp::v_GetTraceBasisKey ( const int  i,
const int  k 
) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 948 of file StdPrismExp.cpp.

950 {
951  ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
952  ASSERTL2(k >= 0 && k <= 1, "basis key id is out of range");
953 
954  switch (i)
955  {
956  case 0:
957  {
959  m_base[k]->GetNumPoints(),
960  m_base[k]->GetNumModes());
961  }
962  case 2:
963  case 4:
964  {
965  return EvaluateQuadFaceBasisKey(k, m_base[k + 1]->GetBasisType(),
966  m_base[k + 1]->GetNumPoints(),
967  m_base[k + 1]->GetNumModes());
968  }
969  case 1:
970  case 3:
971  {
972  return EvaluateTriFaceBasisKey(k, m_base[2 * k]->GetBasisType(),
973  m_base[2 * k]->GetNumPoints(),
974  m_base[2 * k]->GetNumModes());
975  }
976  break;
977  }
978 
979  // Should never get here.
981 }
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::StdPrismExp::v_GetTraceCoeffMap ( const unsigned int  fid,
Array< OneD, unsigned int > &  maparray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1168 of file StdPrismExp.cpp.

1170 {
1172  "Method only implemented if BasisType is identical"
1173  "in x and y directions");
1176  "Method only implemented for Modified_A BasisType"
1177  "(x and y direction) and Modified_B BasisType (z "
1178  "direction)");
1179  int p, q, r, idx = 0;
1180  int P = 0, Q = 0;
1181 
1182  switch (fid)
1183  {
1184  case 0:
1185  P = m_base[0]->GetNumModes();
1186  Q = m_base[1]->GetNumModes();
1187  break;
1188  case 1:
1189  case 3:
1190  P = m_base[0]->GetNumModes();
1191  Q = m_base[2]->GetNumModes();
1192  break;
1193  case 2:
1194  case 4:
1195  P = m_base[1]->GetNumModes();
1196  Q = m_base[2]->GetNumModes();
1197  break;
1198  default:
1199  ASSERTL0(false, "fid must be between 0 and 4");
1200  }
1201 
1202  if (maparray.size() != P * Q)
1203  {
1204  maparray = Array<OneD, unsigned int>(P * Q);
1205  }
1206 
1207  // Set up ordering inside each 2D face. Also for triangular faces,
1208  // populate signarray.
1209  switch (fid)
1210  {
1211  case 0: // Bottom quad
1212  for (q = 0; q < Q; ++q)
1213  {
1214  for (p = 0; p < P; ++p)
1215  {
1216  maparray[q * P + p] = GetMode(p, q, 0);
1217  }
1218  }
1219  break;
1220  case 1: // Left triangle
1221  for (p = 0; p < P; ++p)
1222  {
1223  for (r = 0; r < Q - p; ++r)
1224  {
1225  maparray[idx++] = GetMode(p, 0, r);
1226  }
1227  }
1228  break;
1229  case 2: // Slanted quad
1230  for (q = 0; q < P; ++q)
1231  {
1232  maparray[q] = GetMode(1, q, 0);
1233  }
1234  for (q = 0; q < P; ++q)
1235  {
1236  maparray[P + q] = GetMode(0, q, 1);
1237  }
1238  for (r = 1; r < Q - 1; ++r)
1239  {
1240  for (q = 0; q < P; ++q)
1241  {
1242  maparray[(r + 1) * P + q] = GetMode(1, q, r);
1243  }
1244  }
1245  break;
1246  case 3: // Right triangle
1247  for (p = 0; p < P; ++p)
1248  {
1249  for (r = 0; r < Q - p; ++r)
1250  {
1251  maparray[idx++] = GetMode(p, 1, r);
1252  }
1253  }
1254  break;
1255  case 4: // Rear quad
1256  for (r = 0; r < Q; ++r)
1257  {
1258  for (q = 0; q < P; ++q)
1259  {
1260  maparray[r * P + q] = GetMode(0, q, r);
1261  }
1262  }
1263  break;
1264  default:
1265  ASSERTL0(false, "Face to element map unavailable.");
1266  }
1267 }

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

◆ v_GetTraceInteriorToElementMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1657 of file StdPrismExp.cpp.

1660 {
1661  const int P = m_base[0]->GetNumModes() - 1;
1662  const int Q = m_base[1]->GetNumModes() - 1;
1663  const int R = m_base[2]->GetNumModes() - 1;
1664  const int nFaceIntCoeffs = v_GetTraceIntNcoeffs(fid);
1665  int p, q, r, idx = 0;
1666  int nummodesA = 0;
1667  int nummodesB = 0;
1668  int i = 0;
1669  int j = 0;
1670 
1671  if (maparray.size() != nFaceIntCoeffs)
1672  {
1673  maparray = Array<OneD, unsigned int>(nFaceIntCoeffs);
1674  }
1675 
1676  if (signarray.size() != nFaceIntCoeffs)
1677  {
1678  signarray = Array<OneD, int>(nFaceIntCoeffs, 1);
1679  }
1680  else
1681  {
1682  fill(signarray.get(), signarray.get() + nFaceIntCoeffs, 1);
1683  }
1684 
1685  // Set up an array indexing for quad faces, since the ordering may
1686  // need to be transposed depending on orientation.
1687  Array<OneD, int> arrayindx(nFaceIntCoeffs);
1688  if (fid != 1 && fid != 3)
1689  {
1690  if (fid == 0) // Base quad
1691  {
1692  nummodesA = P - 1;
1693  nummodesB = Q - 1;
1694  }
1695  else // front and back quad
1696  {
1697  nummodesA = Q - 1;
1698  nummodesB = R - 1;
1699  }
1700 
1701  for (i = 0; i < nummodesB; i++)
1702  {
1703  for (j = 0; j < nummodesA; j++)
1704  {
1705  if (faceOrient < eDir1FwdDir2_Dir2FwdDir1)
1706  {
1707  arrayindx[i * nummodesA + j] = i * nummodesA + j;
1708  }
1709  else
1710  {
1711  arrayindx[i * nummodesA + j] = j * nummodesB + i;
1712  }
1713  }
1714  }
1715  }
1716 
1717  switch (fid)
1718  {
1719  case 0: // Bottom quad
1720  for (q = 2; q <= Q; ++q)
1721  {
1722  for (p = 2; p <= P; ++p)
1723  {
1724  maparray[arrayindx[(q - 2) * nummodesA + (p - 2)]] =
1725  GetMode(p, q, 0);
1726  }
1727  }
1728  break;
1729 
1730  case 1: // Left triangle
1731  for (p = 2; p <= P; ++p)
1732  {
1733  for (r = 1; r <= R - p; ++r)
1734  {
1735  if (faceOrient == eDir1BwdDir1_Dir2FwdDir2)
1736  {
1737  signarray[idx] = p % 2 ? -1 : 1;
1738  }
1739  maparray[idx++] = GetMode(p, 0, r);
1740  }
1741  }
1742  break;
1743 
1744  case 2: // Slanted quad
1745  for (r = 1; r <= R - 1; ++r)
1746  {
1747  for (q = 2; q <= Q; ++q)
1748  {
1749  maparray[arrayindx[(r - 1) * nummodesA + (q - 2)]] =
1750  GetMode(1, q, r);
1751  }
1752  }
1753  break;
1754 
1755  case 3: // Right triangle
1756  for (p = 2; p <= P; ++p)
1757  {
1758  for (r = 1; r <= R - p; ++r)
1759  {
1760  if (faceOrient == eDir1BwdDir1_Dir2FwdDir2)
1761  {
1762  signarray[idx] = p % 2 ? -1 : 1;
1763  }
1764  maparray[idx++] = GetMode(p, 1, r);
1765  }
1766  }
1767  break;
1768 
1769  case 4: // Back quad
1770  for (r = 2; r <= R; ++r)
1771  {
1772  for (q = 2; q <= Q; ++q)
1773  {
1774  maparray[arrayindx[(r - 2) * nummodesA + (q - 2)]] =
1775  GetMode(0, q, r);
1776  }
1777  }
1778  break;
1779 
1780  default:
1781  ASSERTL0(false, "Face interior map not available.");
1782  }
1783 
1784  // Triangular faces are processed in the above switch loop; for
1785  // remaining quad faces, set up orientation if necessary.
1786  if (fid == 1 || fid == 3)
1787  return;
1788 
1789  if (faceOrient == eDir1FwdDir1_Dir2BwdDir2 ||
1790  faceOrient == eDir1BwdDir1_Dir2BwdDir2 ||
1791  faceOrient == eDir1BwdDir2_Dir2FwdDir1 ||
1792  faceOrient == eDir1BwdDir2_Dir2BwdDir1)
1793  {
1794  if (faceOrient < eDir1FwdDir2_Dir2FwdDir1)
1795  {
1796  for (i = 1; i < nummodesB; i += 2)
1797  {
1798  for (j = 0; j < nummodesA; j++)
1799  {
1800  signarray[arrayindx[i * nummodesA + j]] *= -1;
1801  }
1802  }
1803  }
1804  else
1805  {
1806  for (i = 0; i < nummodesB; i++)
1807  {
1808  for (j = 1; j < nummodesA; j += 2)
1809  {
1810  signarray[arrayindx[i * nummodesA + j]] *= -1;
1811  }
1812  }
1813  }
1814  }
1815 
1816  if (faceOrient == eDir1BwdDir1_Dir2FwdDir2 ||
1817  faceOrient == eDir1BwdDir1_Dir2BwdDir2 ||
1818  faceOrient == eDir1FwdDir2_Dir2BwdDir1 ||
1819  faceOrient == eDir1BwdDir2_Dir2BwdDir1)
1820  {
1821  if (faceOrient < eDir1FwdDir2_Dir2FwdDir1)
1822  {
1823  for (i = 0; i < nummodesB; i++)
1824  {
1825  for (j = 1; j < nummodesA; j += 2)
1826  {
1827  signarray[arrayindx[i * nummodesA + j]] *= -1;
1828  }
1829  }
1830  }
1831  else
1832  {
1833  for (i = 1; i < nummodesB; i += 2)
1834  {
1835  for (j = 0; j < nummodesA; j++)
1836  {
1837  signarray[arrayindx[i * nummodesA + j]] *= -1;
1838  }
1839  }
1840  }
1841  }
1842 }
virtual int v_GetTraceIntNcoeffs(const int i) const

References ASSERTL0, 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, GetMode(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LibUtilities::P, CellMLToNektar.cellml_metadata::p, and v_GetTraceIntNcoeffs().

◆ v_GetTraceIntNcoeffs()

int Nektar::StdRegions::StdPrismExp::v_GetTraceIntNcoeffs ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 879 of file StdPrismExp.cpp.

880 {
881  ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
882 
883  int Pi = GetBasisNumModes(0) - 2;
884  int Qi = GetBasisNumModes(1) - 2;
885  int Ri = GetBasisNumModes(2) - 2;
886 
887  if (i == 0)
888  {
889  return Pi * Qi;
890  }
891  else if (i == 1 || i == 3)
892  {
893  return Pi * (2 * Ri - Pi - 1) / 2;
894  }
895  else
896  {
897  return Qi * Ri;
898  }
899 }

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

Referenced by v_GetTraceInteriorToElementMap().

◆ v_GetTraceNcoeffs()

int Nektar::StdRegions::StdPrismExp::v_GetTraceNcoeffs ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 861 of file StdPrismExp.cpp.

862 {
863  ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
864  if (i == 0)
865  {
866  return GetBasisNumModes(0) * GetBasisNumModes(1);
867  }
868  else if (i == 1 || i == 3)
869  {
870  int P = GetBasisNumModes(0) - 1, Q = GetBasisNumModes(2) - 1;
871  return Q + 1 + (P * (1 + 2 * Q - P)) / 2;
872  }
873  else
874  {
875  return GetBasisNumModes(1) * GetBasisNumModes(2);
876  }
877 }

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

◆ v_GetTraceNumModes()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 737 of file StdPrismExp.cpp.

739 {
740  int nummodes[3] = {m_base[0]->GetNumModes(), m_base[1]->GetNumModes(),
741  m_base[2]->GetNumModes()};
742  switch (fid)
743  {
744  // base quad
745  case 0:
746  {
747  numModes0 = nummodes[0];
748  numModes1 = nummodes[1];
749  }
750  break;
751  // front and back quad
752  case 2:
753  case 4:
754  {
755  numModes0 = nummodes[1];
756  numModes1 = nummodes[2];
757  }
758  break;
759  // triangles
760  case 1:
761  case 3:
762  {
763  numModes0 = nummodes[0];
764  numModes1 = nummodes[2];
765  }
766  break;
767  }
768 
769  if (faceOrient >= eDir1FwdDir2_Dir2FwdDir1)
770  {
771  std::swap(numModes0, numModes1);
772  }
773 }

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

◆ v_GetTraceNumPoints()

int Nektar::StdRegions::StdPrismExp::v_GetTraceNumPoints ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 910 of file StdPrismExp.cpp.

911 {
912  ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
913 
914  if (i == 0)
915  {
916  return m_base[0]->GetNumPoints() * m_base[1]->GetNumPoints();
917  }
918  else if (i == 1 || i == 3)
919  {
920  return m_base[0]->GetNumPoints() * m_base[2]->GetNumPoints();
921  }
922  else
923  {
924  return m_base[1]->GetNumPoints() * m_base[2]->GetNumPoints();
925  }
926 }

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

◆ v_GetTracePointsKey()

LibUtilities::PointsKey Nektar::StdRegions::StdPrismExp::v_GetTracePointsKey ( const int  i,
const int  j 
) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 928 of file StdPrismExp.cpp.

930 {
931  ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
932  ASSERTL2(j == 0 || j == 1, "face direction is out of range");
933 
934  if (i == 0)
935  {
936  return m_base[j]->GetPointsKey();
937  }
938  else if (i == 1 || i == 3)
939  {
940  return m_base[2 * j]->GetPointsKey();
941  }
942  else
943  {
944  return m_base[j + 1]->GetPointsKey();
945  }
946 }

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

◆ v_GetVertexMap()

int Nektar::StdRegions::StdPrismExp::v_GetVertexMap ( int  localVertexId,
bool  useCoeffPacking = false 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalPrismExp.

Definition at line 1005 of file StdPrismExp.cpp.

1006 {
1010  "Mapping not defined for this type of basis");
1011 
1012  int l = 0;
1013 
1014  if (useCoeffPacking == true) // follow packing of coefficients i.e q,r,p
1015  {
1016  switch (vId)
1017  {
1018  case 0:
1019  l = GetMode(0, 0, 0);
1020  break;
1021  case 1:
1022  l = GetMode(0, 0, 1);
1023  break;
1024  case 2:
1025  l = GetMode(0, 1, 0);
1026  break;
1027  case 3:
1028  l = GetMode(0, 1, 1);
1029  break;
1030  case 4:
1031  l = GetMode(1, 0, 0);
1032  break;
1033  case 5:
1034  l = GetMode(1, 1, 0);
1035  break;
1036  default:
1037  ASSERTL0(false, "local vertex id must be between 0 and 5");
1038  }
1039  }
1040  else
1041  {
1042  switch (vId)
1043  {
1044  case 0:
1045  l = GetMode(0, 0, 0);
1046  break;
1047  case 1:
1048  l = GetMode(1, 0, 0);
1049  break;
1050  case 2:
1051  l = GetMode(1, 1, 0);
1052  break;
1053  case 3:
1054  l = GetMode(0, 1, 0);
1055  break;
1056  case 4:
1057  l = GetMode(0, 0, 1);
1058  break;
1059  case 5:
1060  l = GetMode(0, 1, 1);
1061  break;
1062  default:
1063  ASSERTL0(false, "local vertex id must be between 0 and 5");
1064  }
1065  }
1066 
1067  return l;
1068 }

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

◆ v_IProductWRTBase()

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

Calculate the inner product of inarray with respect to the basis B=base0*base1*base2 and put into outarray:

\( \begin{array}{rcl} I_{pqr} = (\phi_{pqr}, u)_{\delta} & = & \sum_{i=0}^{nq_0} \sum_{j=0}^{nq_1} \sum_{k=0}^{nq_2} \psi_{p}^{a} (\bar \eta_{1i}) \psi_{q}^{a} (\xi_{2j}) \psi_{pr}^{b} (\xi_{3k}) w_i w_j w_k u(\bar \eta_{1,i} \xi_{2,j} \xi_{3,k}) J_{i,j,k}\\ & = & \sum_{i=0}^{nq_0} \psi_p^a(\bar \eta_{1,i}) \sum_{j=0}^{nq_1} \psi_{q}^a(\xi_{2,j}) \sum_{k=0}^{nq_2} \psi_{pr}^b u(\bar \eta_{1i},\xi_{2j},\xi_{3k}) J_{i,j,k} \end{array} \)
where

\( \phi_{pqr} (\xi_1 , \xi_2 , \xi_3) = \psi_p^a (\bar \eta_1) \psi_{q}^a (\xi_2) \psi_{pr}^b (\xi_3) \)
which can be implemented as
\(f_{pr} (\xi_{3k}) = \sum_{k=0}^{nq_3} \psi_{pr}^b u(\bar \eta_{1i},\xi_{2j},\xi_{3k}) J_{i,j,k} = {\bf B_3 U} \)
\( g_{q} (\xi_{3k}) = \sum_{j=0}^{nq_1} \psi_{q}^a (\xi_{2j}) f_{pr} (\xi_{3k}) = {\bf B_2 F} \)
\( (\phi_{pqr}, u)_{\delta} = \sum_{k=0}^{nq_0} \psi_{p}^a (\xi_{3k}) g_{q} (\xi_{3k}) = {\bf B_1 G} \)

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalPrismExp, and Nektar::LocalRegions::PrismExp.

Definition at line 387 of file StdPrismExp.cpp.

389 {
392  "Basis[1] is not a general tensor type");
393 
396  "Basis[2] is not a general tensor type");
397 
398  if (m_base[0]->Collocation() && m_base[1]->Collocation())
399  {
400  MultiplyByQuadratureMetric(inarray, outarray);
401  }
402  else
403  {
404  StdPrismExp::v_IProductWRTBase_SumFac(inarray, outarray);
405  }
406 }
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:731
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)

References ASSERTL1, Nektar::LibUtilities::eModified_B, Nektar::LibUtilities::eModified_C, Nektar::LibUtilities::eOrtho_B, Nektar::LibUtilities::eOrtho_C, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::MultiplyByQuadratureMetric(), and v_IProductWRTBase_SumFac().

Referenced by v_FwdTrans().

◆ v_IProductWRTBase_MatOp()

void Nektar::StdRegions::StdPrismExp::v_IProductWRTBase_MatOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Implementation of the local matrix inner product operation.

Definition at line 411 of file StdPrismExp.cpp.

414 {
415  int nq = GetTotPoints();
416  StdMatrixKey iprodmatkey(eIProductWRTBase, DetShapeType(), *this);
417  DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
418 
419  Blas::Dgemv('N', m_ncoeffs, nq, 1.0, iprodmat->GetPtr().get(), m_ncoeffs,
420  inarray.get(), 1, 0.0, outarray.get(), 1);
421 }
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140
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 = A x where A[m x n].
Definition: Blas.hpp:246

References Nektar::StdRegions::StdExpansion::DetShapeType(), Blas::Dgemv(), Nektar::StdRegions::eIProductWRTBase, Nektar::StdRegions::StdExpansion::GetStdMatrix(), Nektar::StdRegions::StdExpansion::GetTotPoints(), and Nektar::StdRegions::StdExpansion::m_ncoeffs.

◆ v_IProductWRTBase_SumFac()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PrismExp, and Nektar::StdRegions::StdNodalPrismExp.

Definition at line 423 of file StdPrismExp.cpp.

426 {
427  int nquad1 = m_base[1]->GetNumPoints();
428  int nquad2 = m_base[2]->GetNumPoints();
429  int order0 = m_base[0]->GetNumModes();
430  int order1 = m_base[1]->GetNumModes();
431 
432  Array<OneD, NekDouble> wsp(order0 * nquad2 * (nquad1 + order1));
433 
434  if (multiplybyweights)
435  {
436  Array<OneD, NekDouble> tmp(inarray.size());
437 
438  MultiplyByQuadratureMetric(inarray, tmp);
440  m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
441  tmp, outarray, wsp, true, true, true);
442  }
443  else
444  {
446  m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
447  inarray, outarray, wsp, true, true, true);
448  }
449 }
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)

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

Referenced by v_IProductWRTBase(), and Nektar::StdRegions::StdNodalPrismExp::v_IProductWRTBase_SumFac().

◆ v_IProductWRTBase_SumFacKernel()

void Nektar::StdRegions::StdPrismExp::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 
)
protectedvirtual

Implements Nektar::StdRegions::StdExpansion3D.

Definition at line 451 of file StdPrismExp.cpp.

458 {
459  boost::ignore_unused(doCheckCollDir0, doCheckCollDir1, doCheckCollDir2);
460 
461  // Interior prism implementation based on Spen's book page
462  // 119. and 608.
463  const int nquad0 = m_base[0]->GetNumPoints();
464  const int nquad1 = m_base[1]->GetNumPoints();
465  const int nquad2 = m_base[2]->GetNumPoints();
466  const int order0 = m_base[0]->GetNumModes();
467  const int order1 = m_base[1]->GetNumModes();
468  const int order2 = m_base[2]->GetNumModes();
469 
470  int i, mode;
471 
472  ASSERTL1(wsp.size() >= nquad1 * nquad2 * order0 + nquad2 * order0 * order1,
473  "Insufficient workspace size");
474 
475  Array<OneD, NekDouble> tmp0 = wsp;
476  Array<OneD, NekDouble> tmp1 = wsp + nquad1 * nquad2 * order0;
477 
478  // Inner product with respect to the '0' direction
479  Blas::Dgemm('T', 'N', nquad1 * nquad2, order0, nquad0, 1.0, inarray.get(),
480  nquad0, base0.get(), nquad0, 0.0, tmp0.get(), nquad1 * nquad2);
481 
482  // Inner product with respect to the '1' direction
483  Blas::Dgemm('T', 'N', nquad2 * order0, order1, nquad1, 1.0, tmp0.get(),
484  nquad1, base1.get(), nquad1, 0.0, tmp1.get(), nquad2 * order0);
485 
486  // Inner product with respect to the '2' direction
487  for (mode = i = 0; i < order0; ++i)
488  {
489  Blas::Dgemm('T', 'N', order2 - i, order1, nquad2, 1.0,
490  base2.get() + mode * nquad2, nquad2,
491  tmp1.get() + i * nquad2, nquad2 * order0, 0.0,
492  outarray.get() + mode * order1, order2 - i);
493  mode += order2 - i;
494  }
495 
496  // Fix top singular vertices; performs phi_{0,q,1} +=
497  // phi_1(xi_1)*phi_q(xi_2)*phi_{01}*phi_r(xi_2).
499  {
500  for (i = 0; i < order1; ++i)
501  {
502  mode = GetMode(0, i, 1);
503  outarray[mode] +=
504  Blas::Ddot(nquad2, base2.get() + nquad2, 1,
505  tmp1.get() + i * order0 * nquad2 + nquad2, 1);
506  }
507  }
508 }
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:182

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

◆ v_IProductWRTDerivBase()

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

Inner product of inarray over region with respect to the object's default expansion basis; output in outarray.

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalPrismExp, and Nektar::LocalRegions::PrismExp.

Definition at line 514 of file StdPrismExp.cpp.

517 {
518  v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
519 }
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)

References v_IProductWRTDerivBase_SumFac().

◆ v_IProductWRTDerivBase_MatOp()

void Nektar::StdRegions::StdPrismExp::v_IProductWRTDerivBase_MatOp ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Definition at line 521 of file StdPrismExp.cpp.

524 {
525  ASSERTL0(dir >= 0 && dir <= 2, "input dir is out of range");
526 
527  int nq = GetTotPoints();
529 
530  switch (dir)
531  {
532  case 0:
533  mtype = eIProductWRTDerivBase0;
534  break;
535  case 1:
536  mtype = eIProductWRTDerivBase1;
537  break;
538  case 2:
539  mtype = eIProductWRTDerivBase2;
540  break;
541  }
542 
543  StdMatrixKey iprodmatkey(mtype, DetShapeType(), *this);
544  DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
545 
546  Blas::Dgemv('N', m_ncoeffs, nq, 1.0, iprodmat->GetPtr().get(), m_ncoeffs,
547  inarray.get(), 1, 0.0, outarray.get(), 1);
548 }

References ASSERTL0, Nektar::StdRegions::StdExpansion::DetShapeType(), Blas::Dgemv(), Nektar::StdRegions::eIProductWRTDerivBase0, Nektar::StdRegions::eIProductWRTDerivBase1, Nektar::StdRegions::eIProductWRTDerivBase2, Nektar::StdRegions::StdExpansion::GetStdMatrix(), Nektar::StdRegions::StdExpansion::GetTotPoints(), and Nektar::StdRegions::StdExpansion::m_ncoeffs.

◆ v_IProductWRTDerivBase_SumFac()

void Nektar::StdRegions::StdPrismExp::v_IProductWRTDerivBase_SumFac ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalPrismExp, and Nektar::LocalRegions::PrismExp.

Definition at line 550 of file StdPrismExp.cpp.

553 {
554  ASSERTL0(dir >= 0 && dir <= 2, "input dir is out of range");
555 
556  int i;
557  int order0 = m_base[0]->GetNumModes();
558  int order1 = m_base[1]->GetNumModes();
559  int nquad0 = m_base[0]->GetNumPoints();
560  int nquad1 = m_base[1]->GetNumPoints();
561  int nquad2 = m_base[2]->GetNumPoints();
562 
563  const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
564  const Array<OneD, const NekDouble> &z2 = m_base[2]->GetZ();
565  Array<OneD, NekDouble> gfac0(nquad0);
566  Array<OneD, NekDouble> gfac2(nquad2);
567  Array<OneD, NekDouble> tmp0(nquad0 * nquad1 * nquad2);
568  Array<OneD, NekDouble> wsp(order0 * nquad2 * (nquad1 + order1));
569 
570  // set up geometric factor: (1+z0)/2
571  for (i = 0; i < nquad0; ++i)
572  {
573  gfac0[i] = 0.5 * (1 + z0[i]);
574  }
575 
576  // Set up geometric factor: 2/(1-z2)
577  for (i = 0; i < nquad2; ++i)
578  {
579  gfac2[i] = 2.0 / (1 - z2[i]);
580  }
581 
582  // Scale first derivative term by gfac2.
583  if (dir != 1)
584  {
585  for (i = 0; i < nquad2; ++i)
586  {
587  Vmath::Smul(nquad0 * nquad1, gfac2[i],
588  &inarray[0] + i * nquad0 * nquad1, 1,
589  &tmp0[0] + i * nquad0 * nquad1, 1);
590  }
591  MultiplyByQuadratureMetric(tmp0, tmp0);
592  }
593 
594  switch (dir)
595  {
596  case 0:
597  {
599  m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
600  m_base[2]->GetBdata(), tmp0, outarray, wsp, true, true, true);
601  break;
602  }
603 
604  case 1:
605  {
606  MultiplyByQuadratureMetric(inarray, tmp0);
608  m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
609  m_base[2]->GetBdata(), tmp0, outarray, wsp, true, true, true);
610  break;
611  }
612 
613  case 2:
614  {
615  Array<OneD, NekDouble> tmp1(m_ncoeffs);
616 
617  // Scale eta_1 derivative with gfac0.
618  for (i = 0; i < nquad1 * nquad2; ++i)
619  {
620  Vmath::Vmul(nquad0, &gfac0[0], 1, &tmp0[0] + i * nquad0, 1,
621  &tmp0[0] + i * nquad0, 1);
622  }
623 
625  m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
626  m_base[2]->GetBdata(), tmp0, tmp1, wsp, true, true, true);
627 
628  MultiplyByQuadratureMetric(inarray, tmp0);
630  m_base[0]->GetBdata(), m_base[1]->GetBdata(),
631  m_base[2]->GetDbdata(), tmp0, outarray, wsp, true, true, true);
632 
633  Vmath::Vadd(m_ncoeffs, &tmp1[0], 1, &outarray[0], 1, &outarray[0],
634  1);
635  break;
636  }
637  }
638 }
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.cpp:209
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.cpp:359

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

Referenced by v_IProductWRTDerivBase(), and Nektar::StdRegions::StdNodalPrismExp::v_IProductWRTDerivBase_SumFac().

◆ v_IsBoundaryInteriorExpansion()

bool Nektar::StdRegions::StdPrismExp::v_IsBoundaryInteriorExpansion ( )
protectedvirtual

◆ v_LocCollapsedToLocCoord()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 664 of file StdPrismExp.cpp.

666 {
667  xi[0] = (1.0 + eta[0]) * (1.0 - eta[2]) * 0.5 - 1.0;
668  xi[1] = eta[1];
669  xi[2] = eta[2];
670 }

◆ v_LocCoordToLocCollapsed()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 644 of file StdPrismExp.cpp.

646 {
647  NekDouble d2 = 1.0 - xi[2];
648  if (fabs(d2) < NekConstants::kNekZeroTol)
649  {
650  if (d2 >= 0.)
651  {
653  }
654  else
655  {
657  }
658  }
659  eta[2] = xi[2]; // eta_z = xi_z
660  eta[1] = xi[1]; // eta_y = xi_y
661  eta[0] = 2.0 * (1.0 + xi[0]) / d2 - 1.0;
662 }
static const NekDouble kNekZeroTol

References Nektar::NekConstants::kNekZeroTol.

◆ v_MultiplyByStdQuadratureMetric()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1969 of file StdPrismExp.cpp.

1972 {
1973  int i, j;
1974  int nquad0 = m_base[0]->GetNumPoints();
1975  int nquad1 = m_base[1]->GetNumPoints();
1976  int nquad2 = m_base[2]->GetNumPoints();
1977 
1978  const Array<OneD, const NekDouble> &w0 = m_base[0]->GetW();
1979  const Array<OneD, const NekDouble> &w1 = m_base[1]->GetW();
1980  const Array<OneD, const NekDouble> &w2 = m_base[2]->GetW();
1981 
1982  const Array<OneD, const NekDouble> &z2 = m_base[2]->GetZ();
1983 
1984  // Multiply by integration constants in x-direction
1985  for (i = 0; i < nquad1 * nquad2; ++i)
1986  {
1987  Vmath::Vmul(nquad0, inarray.get() + i * nquad0, 1, w0.get(), 1,
1988  outarray.get() + i * nquad0, 1);
1989  }
1990 
1991  // Multiply by integration constants in y-direction
1992  for (j = 0; j < nquad2; ++j)
1993  {
1994  for (i = 0; i < nquad1; ++i)
1995  {
1996  Blas::Dscal(nquad0, w1[i],
1997  &outarray[0] + i * nquad0 + j * nquad0 * nquad1, 1);
1998  }
1999  }
2000 
2001  // Multiply by integration constants in z-direction; need to
2002  // incorporate factor (1-eta_3)/2 into weights, but only if using
2003  // GLL quadrature points.
2004  switch (m_base[2]->GetPointsType())
2005  {
2006  // (1,0) Jacobi inner product.
2007  case LibUtilities::eGaussRadauMAlpha1Beta0:
2008  for (i = 0; i < nquad2; ++i)
2009  {
2010  Blas::Dscal(nquad0 * nquad1, 0.5 * w2[i],
2011  &outarray[0] + i * nquad0 * nquad1, 1);
2012  }
2013  break;
2014 
2015  default:
2016  for (i = 0; i < nquad2; ++i)
2017  {
2018  Blas::Dscal(nquad0 * nquad1, 0.5 * (1 - z2[i]) * w2[i],
2019  &outarray[0] + i * nquad0 * nquad1, 1);
2020  }
2021  break;
2022  }
2023 }
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:213
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition: Blas.hpp:168

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

◆ v_NumBndryCoeffs()

int Nektar::StdRegions::StdPrismExp::v_NumBndryCoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 821 of file StdPrismExp.cpp.

822 {
825  "BasisType is not a boundary interior form");
828  "BasisType is not a boundary interior form");
831  "BasisType is not a boundary interior form");
832 
833  int P = m_base[0]->GetNumModes();
834  int Q = m_base[1]->GetNumModes();
835  int R = m_base[2]->GetNumModes();
836 
838 }
int getNumberOfBndCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:293

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

◆ v_NumDGBndryCoeffs()

int Nektar::StdRegions::StdPrismExp::v_NumDGBndryCoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 840 of file StdPrismExp.cpp.

841 {
844  "BasisType is not a boundary interior form");
847  "BasisType is not a boundary interior form");
850  "BasisType is not a boundary interior form");
851 
852  int P = m_base[0]->GetNumModes() - 1;
853  int Q = m_base[1]->GetNumModes() - 1;
854  int R = m_base[2]->GetNumModes() - 1;
855 
856  return (P + 1) * (Q + 1) // 1 rect. face on base
857  + 2 * (Q + 1) * (R + 1) // other 2 rect. faces
858  + 2 * (R + 1) + P * (1 + 2 * R - P); // 2 tri. faces
859 }

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

◆ v_PhysDeriv() [1/2]

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

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::PrismExp.

Definition at line 94 of file StdPrismExp.cpp.

98 {
99  int Qx = m_base[0]->GetNumPoints();
100  int Qy = m_base[1]->GetNumPoints();
101  int Qz = m_base[2]->GetNumPoints();
102  int Qtot = Qx * Qy * Qz;
103 
104  Array<OneD, NekDouble> dEta_bar1(Qtot, 0.0);
105 
106  Array<OneD, const NekDouble> eta_x, eta_z;
107  eta_x = m_base[0]->GetZ();
108  eta_z = m_base[2]->GetZ();
109 
110  int i, k;
111 
112  bool Do_1 = (out_dxi1.size() > 0) ? true : false;
113  bool Do_3 = (out_dxi3.size() > 0) ? true : false;
114 
115  // out_dXi2 is just a tensor derivative so is just passed through
116  if (Do_3)
117  {
118  PhysTensorDeriv(u_physical, dEta_bar1, out_dxi2, out_dxi3);
119  }
120  else if (Do_1)
121  {
122  PhysTensorDeriv(u_physical, dEta_bar1, out_dxi2, NullNekDouble1DArray);
123  }
124  else // case if just require 2nd direction
125  {
126  PhysTensorDeriv(u_physical, NullNekDouble1DArray, out_dxi2,
128  }
129 
130  if (Do_1)
131  {
132  for (k = 0; k < Qz; ++k)
133  {
134  Vmath::Smul(Qx * Qy, 2.0 / (1.0 - eta_z[k]),
135  &dEta_bar1[0] + k * Qx * Qy, 1,
136  &out_dxi1[0] + k * Qx * Qy, 1);
137  }
138  }
139 
140  if (Do_3)
141  {
142  // divide dEta_Bar1 by (1-eta_z)
143  for (k = 0; k < Qz; ++k)
144  {
145  Vmath::Smul(Qx * Qy, 1.0 / (1.0 - eta_z[k]),
146  &dEta_bar1[0] + k * Qx * Qy, 1,
147  &dEta_bar1[0] + k * Qx * Qy, 1);
148  }
149 
150  // Multiply dEta_Bar1 by (1+eta_x) and add ot out_dxi3
151  for (i = 0; i < Qx; ++i)
152  {
153  Vmath::Svtvp(Qz * Qy, 1.0 + eta_x[i], &dEta_bar1[0] + i, Qx,
154  &out_dxi3[0] + i, Qx, &out_dxi3[0] + i, Qx);
155  }
156  }
157 }
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.
static Array< OneD, NekDouble > NullNekDouble1DArray
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
Definition: Vmath.cpp:622

References Nektar::StdRegions::StdExpansion::m_base, Nektar::NullNekDouble1DArray, Nektar::StdRegions::StdExpansion3D::PhysTensorDeriv(), Vmath::Smul(), and Vmath::Svtvp().

Referenced by v_PhysDeriv(), and v_StdPhysDeriv().

◆ v_PhysDeriv() [2/2]

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

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

See also
StdRegions::StdExpansion::PhysDeriv

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 159 of file StdPrismExp.cpp.

162 {
163  switch (dir)
164  {
165  case 0:
166  {
167  v_PhysDeriv(inarray, outarray, NullNekDouble1DArray,
169  break;
170  }
171 
172  case 1:
173  {
174  v_PhysDeriv(inarray, NullNekDouble1DArray, outarray,
176  break;
177  }
178 
179  case 2:
180  {
182  outarray);
183  break;
184  }
185 
186  default:
187  {
188  ASSERTL1(false, "input dir is out of range");
189  }
190  break;
191  }
192 }
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
Calculate the derivative of the physical points.
Definition: StdPrismExp.cpp:94

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

◆ v_PhysEvaluateBasis()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 706 of file StdPrismExp.cpp.

708 {
709  Array<OneD, NekDouble> coll(3);
710  LocCoordToLocCollapsed(coords, coll);
711 
712  const int nm1 = m_base[1]->GetNumModes();
713  const int nm2 = m_base[2]->GetNumModes();
714  const int b = 2 * nm2 + 1;
715 
716  const int mode0 = floor(0.5 * (b - sqrt(b * b - 8.0 * mode / nm1)));
717  const int tmp =
718  mode - nm1 * (mode0 * (nm2 - 1) + 1 - (mode0 - 2) * (mode0 - 1) / 2);
719  const int mode1 = tmp / (nm2 - mode0);
720  const int mode2 = tmp % (nm2 - mode0);
721 
722  if (mode0 == 0 && mode2 == 1 &&
724  {
725  // handle collapsed top edge to remove mode0 terms
726  return StdExpansion::BaryEvaluateBasis<1>(coll[1], mode1) *
727  StdExpansion::BaryEvaluateBasis<2>(coll[2], mode2);
728  }
729  else
730  {
731  return StdExpansion::BaryEvaluateBasis<0>(coll[0], mode0) *
732  StdExpansion::BaryEvaluateBasis<1>(coll[1], mode1) *
733  StdExpansion::BaryEvaluateBasis<2>(coll[2], mode2);
734  }
735 }
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:291

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

◆ v_ReduceOrderCoeffs()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2163 of file StdPrismExp.cpp.

2166 {
2167  int nquad0 = m_base[0]->GetNumPoints();
2168  int nquad1 = m_base[1]->GetNumPoints();
2169  int nquad2 = m_base[2]->GetNumPoints();
2170  int nqtot = nquad0 * nquad1 * nquad2;
2171  int nmodes0 = m_base[0]->GetNumModes();
2172  int nmodes1 = m_base[1]->GetNumModes();
2173  int nmodes2 = m_base[2]->GetNumModes();
2174  int numMax = nmodes0;
2175 
2176  Array<OneD, NekDouble> coeff(m_ncoeffs);
2177  Array<OneD, NekDouble> coeff_tmp1(m_ncoeffs, 0.0);
2178  Array<OneD, NekDouble> phys_tmp(nqtot, 0.0);
2179  Array<OneD, NekDouble> tmp, tmp2, tmp3, tmp4;
2180 
2181  const LibUtilities::PointsKey Pkey0 = m_base[0]->GetPointsKey();
2182  const LibUtilities::PointsKey Pkey1 = m_base[1]->GetPointsKey();
2183  const LibUtilities::PointsKey Pkey2 = m_base[2]->GetPointsKey();
2184 
2185  LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A, nmodes0, Pkey0);
2186  LibUtilities::BasisKey bortho1(LibUtilities::eOrtho_A, nmodes1, Pkey1);
2187  LibUtilities::BasisKey bortho2(LibUtilities::eOrtho_B, nmodes2, Pkey2);
2188 
2189  int cnt = 0;
2190  int u = 0;
2191  int i = 0;
2192  StdRegions::StdPrismExpSharedPtr OrthoPrismExp;
2193 
2195  bortho0, bortho1, bortho2);
2196 
2197  BwdTrans(inarray, phys_tmp);
2198  OrthoPrismExp->FwdTrans(phys_tmp, coeff);
2199 
2200  // filtering
2201  for (u = 0; u < numMin; ++u)
2202  {
2203  for (i = 0; i < numMin; ++i)
2204  {
2205  Vmath::Vcopy(numMin - u, tmp = coeff + cnt, 1,
2206  tmp2 = coeff_tmp1 + cnt, 1);
2207  cnt += numMax - u;
2208  }
2209 
2210  for (i = numMin; i < numMax; ++i)
2211  {
2212  cnt += numMax - u;
2213  }
2214  }
2215 
2216  OrthoPrismExp->BwdTrans(coeff_tmp1, phys_tmp);
2217  StdPrismExp::FwdTrans(phys_tmp, outarray);
2218 }
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:432
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:44
std::shared_ptr< StdPrismExp > StdPrismExpSharedPtr
Definition: StdPrismExp.h:231

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

◆ v_StdPhysDeriv() [1/2]

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 194 of file StdPrismExp.cpp.

198 {
199  StdPrismExp::v_PhysDeriv(inarray, out_d0, out_d1, out_d2);
200 }

References v_PhysDeriv().

◆ v_StdPhysDeriv() [2/2]

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 202 of file StdPrismExp.cpp.

205 {
206  StdPrismExp::v_PhysDeriv(dir, inarray, outarray);
207 }

References v_PhysDeriv().

◆ v_SVVLaplacianFilter()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PrismExp.

Definition at line 2025 of file StdPrismExp.cpp.

2027 {
2028  // Generate an orthonogal expansion
2029  int qa = m_base[0]->GetNumPoints();
2030  int qb = m_base[1]->GetNumPoints();
2031  int qc = m_base[2]->GetNumPoints();
2032  int nmodes_a = m_base[0]->GetNumModes();
2033  int nmodes_b = m_base[1]->GetNumModes();
2034  int nmodes_c = m_base[2]->GetNumModes();
2035  // Declare orthogonal basis.
2036  LibUtilities::PointsKey pa(qa, m_base[0]->GetPointsType());
2037  LibUtilities::PointsKey pb(qb, m_base[1]->GetPointsType());
2038  LibUtilities::PointsKey pc(qc, m_base[2]->GetPointsType());
2039 
2040  LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A, nmodes_a, pa);
2041  LibUtilities::BasisKey Bb(LibUtilities::eOrtho_A, nmodes_b, pb);
2042  LibUtilities::BasisKey Bc(LibUtilities::eOrtho_B, nmodes_c, pc);
2043  StdPrismExp OrthoExp(Ba, Bb, Bc);
2044 
2045  Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
2046  int i, j, k, cnt = 0;
2047 
2048  // project onto modal space.
2049  OrthoExp.FwdTrans(array, orthocoeffs);
2050 
2051  if (mkey.ConstFactorExists(eFactorSVVPowerKerDiffCoeff))
2052  {
2053  // Rodrigo's power kernel
2054  NekDouble cutoff = mkey.GetConstFactor(eFactorSVVCutoffRatio);
2055  NekDouble SvvDiffCoeff =
2056  mkey.GetConstFactor(eFactorSVVPowerKerDiffCoeff) *
2057  mkey.GetConstFactor(eFactorSVVDiffCoeff);
2058 
2059  for (int i = 0; i < nmodes_a; ++i)
2060  {
2061  for (int j = 0; j < nmodes_b; ++j)
2062  {
2063  NekDouble fac1 = std::max(
2064  pow((1.0 * i) / (nmodes_a - 1), cutoff * nmodes_a),
2065  pow((1.0 * j) / (nmodes_b - 1), cutoff * nmodes_b));
2066 
2067  for (int k = 0; k < nmodes_c - i; ++k)
2068  {
2069  NekDouble fac =
2070  std::max(fac1, pow((1.0 * k) / (nmodes_c - 1),
2071  cutoff * nmodes_c));
2072 
2073  orthocoeffs[cnt] *= SvvDiffCoeff * fac;
2074  cnt++;
2075  }
2076  }
2077  }
2078  }
2079  else if (mkey.ConstFactorExists(
2080  eFactorSVVDGKerDiffCoeff)) // Rodrigo/Mansoor's DG Kernel
2081  {
2082  NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDGKerDiffCoeff) *
2083  mkey.GetConstFactor(eFactorSVVDiffCoeff);
2084 
2085  int max_abc = max(nmodes_a - kSVVDGFiltermodesmin,
2086  nmodes_b - kSVVDGFiltermodesmin);
2087  max_abc = max(max_abc, nmodes_c - kSVVDGFiltermodesmin);
2088  // clamp max_abc
2089  max_abc = max(max_abc, 0);
2090  max_abc = min(max_abc, kSVVDGFiltermodesmax - kSVVDGFiltermodesmin);
2091 
2092  for (int i = 0; i < nmodes_a; ++i)
2093  {
2094  for (int j = 0; j < nmodes_b; ++j)
2095  {
2096  int maxij = max(i, j);
2097 
2098  for (int k = 0; k < nmodes_c - i; ++k)
2099  {
2100  int maxijk = max(maxij, k);
2101  maxijk = min(maxijk, kSVVDGFiltermodesmax - 1);
2102 
2103  orthocoeffs[cnt] *=
2104  SvvDiffCoeff * kSVVDGFilter[max_abc][maxijk];
2105  cnt++;
2106  }
2107  }
2108  }
2109  }
2110  else
2111  {
2112  // SVV filter paramaters (how much added diffusion relative
2113  // to physical one and fraction of modes from which you
2114  // start applying this added diffusion)
2115  //
2116  NekDouble SvvDiffCoeff =
2117  mkey.GetConstFactor(StdRegions::eFactorSVVDiffCoeff);
2118  NekDouble SVVCutOff =
2119  mkey.GetConstFactor(StdRegions::eFactorSVVCutoffRatio);
2120 
2121  // Defining the cut of mode
2122  int cutoff_a = (int)(SVVCutOff * nmodes_a);
2123  int cutoff_b = (int)(SVVCutOff * nmodes_b);
2124  int cutoff_c = (int)(SVVCutOff * nmodes_c);
2125  // To avoid the fac[j] from blowing up
2126  NekDouble epsilon = 1;
2127 
2128  int nmodes = min(min(nmodes_a, nmodes_b), nmodes_c);
2129  NekDouble cutoff = min(min(cutoff_a, cutoff_b), cutoff_c);
2130 
2131  //------"New" Version August 22nd '13--------------------
2132  for (i = 0; i < nmodes_a; ++i) // P
2133  {
2134  for (j = 0; j < nmodes_b; ++j) // Q
2135  {
2136  for (k = 0; k < nmodes_c - i; ++k) // R
2137  {
2138  if (j >= cutoff || i + k >= cutoff)
2139  {
2140  orthocoeffs[cnt] *=
2141  (SvvDiffCoeff *
2142  exp(-(i + k - nmodes) * (i + k - nmodes) /
2143  ((NekDouble)((i + k - cutoff + epsilon) *
2144  (i + k - cutoff + epsilon)))) *
2145  exp(-(j - nmodes) * (j - nmodes) /
2146  ((NekDouble)((j - cutoff + epsilon) *
2147  (j - cutoff + epsilon)))));
2148  }
2149  else
2150  {
2151  orthocoeffs[cnt] *= 0.0;
2152  }
2153  cnt++;
2154  }
2155  }
2156  }
2157  }
2158 
2159  // backward transform to physical space
2160  OrthoExp.BwdTrans(orthocoeffs, array);
2161 }
const int kSVVDGFiltermodesmin
Definition: StdRegions.hpp:352
const int kSVVDGFiltermodesmax
Definition: StdRegions.hpp:353
const NekDouble kSVVDGFilter[9][11]
Definition: StdRegions.hpp:355

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::eOrtho_B, 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.