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

Protected Member Functions

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

Private Member Functions

int GetMode (int I, int J, int K)
 Compute the 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 44 of file StdPrismExp.h.

Constructor & Destructor Documentation

◆ StdPrismExp() [1/3]

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

Definition at line 43 of file StdPrismExp.cpp.

47 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
48 3, Ba, Bb, Bc),
50 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
51 Ba, Bb, Bc)
52{
53 ASSERTL0(Ba.GetNumModes() <= Bc.GetNumModes(),
54 "order in 'a' direction is higher than order in 'c' direction");
55}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
StdExpansion()
Default Constructor.
int getNumberOfCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:279

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

◆ StdPrismExp() [2/3]

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

◆ StdPrismExp() [3/3]

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

◆ ~StdPrismExp()

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

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 1935 of file StdPrismExp.cpp.

1936{
1937 int Q = m_base[1]->GetNumModes() - 1;
1938 int R = m_base[2]->GetNumModes() - 1;
1939
1940 return r + // Skip along stacks (r-direction)
1941 q * (R + 1 - p) + // Skip along columns (q-direction)
1942 (Q + 1) * (p * R + 1 -
1943 (p - 2) * (p - 1) / 2); // Skip along rows (p-direction)
1944}
Array< OneD, LibUtilities::BasisSharedPtr > m_base
std::vector< double > q(NPUPPER *NPUPPER)

References Nektar::StdRegions::StdExpansion::m_base, CellMLToNektar.cellml_metadata::p, and Nektar::UnitTests::q().

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 
)
overrideprotectedvirtual
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.

Definition at line 220 of file StdPrismExp.cpp.

222{
225 "Basis[1] is not a general tensor type");
226
229 "Basis[2] is not a general tensor type");
230
231 if (m_base[0]->Collocation() && m_base[1]->Collocation() &&
232 m_base[2]->Collocation())
233 {
235 m_base[2]->GetNumPoints(),
236 inarray, 1, outarray, 1);
237 }
238 else
239 {
240 StdPrismExp::v_BwdTrans_SumFac(inarray, outarray);
241 }
242}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:156
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:218
void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:49
@ eModified_C
Principle Modified Functions .
Definition: BasisType.h:50
@ eOrtho_C
Principle Orthogonal Functions .
Definition: BasisType.h:46
@ eOrtho_B
Principle Orthogonal Functions .
Definition: BasisType.h:44
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

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 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 244 of file StdPrismExp.cpp.

246{
247 int nquad1 = m_base[1]->GetNumPoints();
248 int nquad2 = m_base[2]->GetNumPoints();
249 int order0 = m_base[0]->GetNumModes();
250 int order1 = m_base[1]->GetNumModes();
251
252 Array<OneD, NekDouble> wsp(nquad2 * order1 * order0 +
253 nquad1 * nquad2 * order0);
254
255 BwdTrans_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetBdata(),
256 m_base[2]->GetBdata(), inarray, outarray, wsp, true,
257 true, true);
258}
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 
)
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion3D.

Definition at line 260 of file StdPrismExp.cpp.

269{
270 int i, mode;
271 int nquad0 = m_base[0]->GetNumPoints();
272 int nquad1 = m_base[1]->GetNumPoints();
273 int nquad2 = m_base[2]->GetNumPoints();
274 int nummodes0 = m_base[0]->GetNumModes();
275 int nummodes1 = m_base[1]->GetNumModes();
276 int nummodes2 = m_base[2]->GetNumModes();
277 Array<OneD, NekDouble> tmp0 = wsp;
278 Array<OneD, NekDouble> tmp1 = tmp0 + nquad2 * nummodes1 * nummodes0;
279
280 for (i = mode = 0; i < nummodes0; ++i)
281 {
282 Blas::Dgemm('N', 'N', nquad2, nummodes1, nummodes2 - i, 1.0,
283 base2.get() + mode * nquad2, nquad2,
284 inarray.get() + mode * nummodes1, nummodes2 - i, 0.0,
285 tmp0.get() + i * nquad2 * nummodes1, nquad2);
286 mode += nummodes2 - i;
287 }
288
290 {
291 for (i = 0; i < nummodes1; i++)
292 {
293 Blas::Daxpy(nquad2, inarray[1 + i * nummodes2],
294 base2.get() + nquad2, 1,
295 tmp0.get() + nquad2 * (nummodes1 + i), 1);
296 }
297 }
298
299 for (i = 0; i < nummodes0; i++)
300 {
301 Blas::Dgemm('N', 'T', nquad1, nquad2, nummodes1, 1.0, base1.get(),
302 nquad1, tmp0.get() + i * nquad2 * nummodes1, nquad2, 0.0,
303 tmp1.get() + i * nquad2 * nquad1, nquad1);
304 }
305
306 Blas::Dgemm('N', 'T', nquad0, nquad2 * nquad1, nummodes0, 1.0, base0.get(),
307 nquad0, tmp1.get(), nquad2 * nquad1, 0.0, outarray.get(),
308 nquad0);
309}
static void Dgemm(const char &transa, const char &transb, const int &m, const int &n, const int &k, const double &alpha, const double *a, const int &lda, const double *b, const int &ldb, const double &beta, double *c, const int &ldc)
BLAS level 3: Matrix-matrix multiply C = A x B where op(A)[m x k], op(B)[k x n], C[m x n] DGEMM perfo...
Definition: Blas.hpp:383
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
Definition: Blas.hpp:135
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48

References Blas::Daxpy(), Blas::Dgemm(), 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 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 958 of file StdPrismExp.cpp.

960{
962 nummodes[modes_offset], nummodes[modes_offset + 1],
963 nummodes[modes_offset + 2]);
964
965 modes_offset += 3;
966 return nmodes;
967}

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

◆ v_CreateStdMatrix()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PrismExp.

Definition at line 1912 of file StdPrismExp.cpp.

1913{
1914 return v_GenMatrix(mkey);
1915}
DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override

References v_GenMatrix().

◆ v_DetShapeType()

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

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 800 of file StdPrismExp.cpp.

801{
803}

References Nektar::LibUtilities::ePrism.

◆ v_FillMode()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 683 of file StdPrismExp.cpp.

684{
685 Array<OneD, NekDouble> tmp(m_ncoeffs, 0.0);
686 tmp[mode] = 1.0;
687 StdPrismExp::v_BwdTrans(tmp, outarray);
688}
void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override

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 
)
overrideprotectedvirtual

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

Inputs:

  • inarray: array of physical quadrature points to be transformed

Outputs:

  • outarray: updated array of expansion coefficients.

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PrismExp.

Definition at line 322 of file StdPrismExp.cpp.

324{
325 v_IProductWRTBase(inarray, outarray);
326
327 // Get Mass matrix inverse
328 StdMatrixKey masskey(eInvMass, DetShapeType(), *this);
329 DNekMatSharedPtr matsys = GetStdMatrix(masskey);
330
331 // copy inarray in case inarray == outarray
332 DNekVec in(m_ncoeffs, outarray);
333 DNekVec out(m_ncoeffs, outarray, eWrapper);
334
335 out = (*matsys) * in;
336}
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:603
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:367
void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
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)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PrismExp.

Definition at line 1825 of file StdPrismExp.cpp.

1826{
1827
1828 MatrixType mtype = mkey.GetMatrixType();
1829
1830 DNekMatSharedPtr Mat;
1831
1832 switch (mtype)
1833 {
1835 {
1836 int nq0 = m_base[0]->GetNumPoints();
1837 int nq1 = m_base[1]->GetNumPoints();
1838 int nq2 = m_base[2]->GetNumPoints();
1839 int nq;
1840
1841 // take definition from key
1842 if (mkey.ConstFactorExists(eFactorConst))
1843 {
1844 nq = (int)mkey.GetConstFactor(eFactorConst);
1845 }
1846 else
1847 {
1848 nq = max(nq0, max(nq1, nq2));
1849 }
1850
1851 int neq =
1853 Array<OneD, Array<OneD, NekDouble>> coords(neq);
1854 Array<OneD, NekDouble> coll(3);
1855 Array<OneD, DNekMatSharedPtr> I(3);
1856 Array<OneD, NekDouble> tmp(nq0);
1857
1858 Mat =
1859 MemoryManager<DNekMat>::AllocateSharedPtr(neq, nq0 * nq1 * nq2);
1860 int cnt = 0;
1861 for (int i = 0; i < nq; ++i)
1862 {
1863 for (int j = 0; j < nq; ++j)
1864 {
1865 for (int k = 0; k < nq - i; ++k, ++cnt)
1866 {
1867 coords[cnt] = Array<OneD, NekDouble>(3);
1868 coords[cnt][0] = -1.0 + 2 * k / (NekDouble)(nq - 1);
1869 coords[cnt][1] = -1.0 + 2 * j / (NekDouble)(nq - 1);
1870 coords[cnt][2] = -1.0 + 2 * i / (NekDouble)(nq - 1);
1871 }
1872 }
1873 }
1874
1875 for (int i = 0; i < neq; ++i)
1876 {
1877 LocCoordToLocCollapsed(coords[i], coll);
1878
1879 I[0] = m_base[0]->GetI(coll);
1880 I[1] = m_base[1]->GetI(coll + 1);
1881 I[2] = m_base[2]->GetI(coll + 2);
1882
1883 // interpolate first coordinate direction
1884 NekDouble fac;
1885 for (int k = 0; k < nq2; ++k)
1886 {
1887 for (int j = 0; j < nq1; ++j)
1888 {
1889
1890 fac = (I[1]->GetPtr())[j] * (I[2]->GetPtr())[k];
1891 Vmath::Smul(nq0, fac, I[0]->GetPtr(), 1, tmp, 1);
1892
1893 Vmath::Vcopy(nq0, &tmp[0], 1,
1894 Mat->GetRawPtr() + k * nq0 * nq1 * neq +
1895 j * nq0 * neq + i,
1896 neq);
1897 }
1898 }
1899 }
1900 }
1901 break;
1902 default:
1903 {
1905 }
1906 break;
1907 }
1908
1909 return Mat;
1910}
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.
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.hpp:100

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)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1083 of file StdPrismExp.cpp.

1084{
1087 "BasisType is not a boundary interior form");
1090 "BasisType is not a boundary interior form");
1093 "BasisType is not a boundary interior form");
1094
1095 int P = m_base[0]->GetNumModes() - 1, p;
1096 int Q = m_base[1]->GetNumModes() - 1, q;
1097 int R = m_base[2]->GetNumModes() - 1, r;
1098 int idx = 0;
1099
1100 int nBnd = NumBndryCoeffs();
1101
1102 if (maparray.size() != nBnd)
1103 {
1104 maparray = Array<OneD, unsigned int>(nBnd);
1105 }
1106
1107 // Loop over all boundary modes (in ascending order).
1108 for (p = 0; p <= P; ++p)
1109 {
1110 // First two q-r planes are entirely boundary modes.
1111 if (p <= 1)
1112 {
1113 for (q = 0; q <= Q; ++q)
1114 {
1115 for (r = 0; r <= R - p; ++r)
1116 {
1117 maparray[idx++] = GetMode(p, q, r);
1118 }
1119 }
1120 }
1121 else
1122 {
1123 // Remaining q-r planes contain boundary modes on the two
1124 // left-hand sides and bottom edge.
1125 for (q = 0; q <= Q; ++q)
1126 {
1127 if (q <= 1)
1128 {
1129 for (r = 0; r <= R - p; ++r)
1130 {
1131 maparray[idx++] = GetMode(p, q, r);
1132 }
1133 }
1134 else
1135 {
1136 maparray[idx++] = GetMode(p, q, 0);
1137 }
1138 }
1139 }
1140 }
1141}
int GetMode(int I, int J, int K)
Compute the local mode number in the expansion for a particular tensorial combination.
@ P
Monomial polynomials .
Definition: BasisType.h:62
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:56

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, CellMLToNektar.cellml_metadata::p, and Nektar::UnitTests::q().

◆ v_GetCoords()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PrismExp.

Definition at line 608 of file StdPrismExp.cpp.

611{
612 Array<OneD, const NekDouble> etaBar_x = m_base[0]->GetZ();
613 Array<OneD, const NekDouble> eta_y = m_base[1]->GetZ();
614 Array<OneD, const NekDouble> eta_z = m_base[2]->GetZ();
615 int Qx = GetNumPoints(0);
616 int Qy = GetNumPoints(1);
617 int Qz = GetNumPoints(2);
618
619 // Convert collapsed coordinates into cartesian coordinates: eta --> xi
620 for (int k = 0; k < Qz; ++k)
621 {
622 for (int j = 0; j < Qy; ++j)
623 {
624 for (int i = 0; i < Qx; ++i)
625 {
626 int s = i + Qx * (j + Qy * k);
627 xi_x[s] = (1.0 - eta_z[k]) * (1.0 + etaBar_x[i]) / 2.0 - 1.0;
628 xi_y[s] = eta_y[j];
629 xi_z[s] = eta_z[k];
630 }
631 }
632 }
633}

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 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 1520 of file StdPrismExp.cpp.

1523{
1524 int i;
1525 bool signChange;
1526 const int P = m_base[0]->GetNumModes() - 1;
1527 const int Q = m_base[1]->GetNumModes() - 1;
1528 const int R = m_base[2]->GetNumModes() - 1;
1529 const int nEdgeIntCoeffs = v_GetEdgeNcoeffs(eid) - 2;
1530
1531 if (maparray.size() != nEdgeIntCoeffs)
1532 {
1533 maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
1534 }
1535
1536 if (signarray.size() != nEdgeIntCoeffs)
1537 {
1538 signarray = Array<OneD, int>(nEdgeIntCoeffs, 1);
1539 }
1540 else
1541 {
1542 fill(signarray.get(), signarray.get() + nEdgeIntCoeffs, 1);
1543 }
1544
1545 // If edge is oriented backwards, change sign of modes which have
1546 // degree 2n+1, n >= 1.
1547 signChange = edgeOrient == eBackwards;
1548
1549 switch (eid)
1550 {
1551 case 0:
1552 for (i = 2; i <= P; ++i)
1553 {
1554 maparray[i - 2] = GetMode(i, 0, 0);
1555 }
1556 break;
1557
1558 case 1:
1559 for (i = 2; i <= Q; ++i)
1560 {
1561 maparray[i - 2] = GetMode(1, i, 0);
1562 }
1563 break;
1564
1565 case 2:
1566 // Base quad; reverse direction.
1567 // signChange = !signChange;
1568 for (i = 2; i <= P; ++i)
1569 {
1570 maparray[i - 2] = GetMode(i, 1, 0);
1571 }
1572 break;
1573
1574 case 3:
1575 // Base quad; reverse direction.
1576 // signChange = !signChange;
1577 for (i = 2; i <= Q; ++i)
1578 {
1579 maparray[i - 2] = GetMode(0, i, 0);
1580 }
1581 break;
1582
1583 case 4:
1584 for (i = 2; i <= R; ++i)
1585 {
1586 maparray[i - 2] = GetMode(0, 0, i);
1587 }
1588 break;
1589
1590 case 5:
1591 for (i = 1; i <= R - 1; ++i)
1592 {
1593 maparray[i - 1] = GetMode(1, 0, i);
1594 }
1595 break;
1596
1597 case 6:
1598 for (i = 1; i <= R - 1; ++i)
1599 {
1600 maparray[i - 1] = GetMode(1, 1, i);
1601 }
1602 break;
1603
1604 case 7:
1605 for (i = 2; i <= R; ++i)
1606 {
1607 maparray[i - 2] = GetMode(0, 1, i);
1608 }
1609 break;
1610
1611 case 8:
1612 for (i = 2; i <= Q; ++i)
1613 {
1614 maparray[i - 2] = GetMode(0, i, 1);
1615 }
1616 break;
1617
1618 default:
1619 ASSERTL0(false, "Edge not defined.");
1620 break;
1621 }
1622
1623 if (signChange)
1624 {
1625 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1626 {
1627 signarray[i] = -1;
1628 }
1629 }
1630}
int v_GetEdgeNcoeffs(const int i) const override

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
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 759 of file StdPrismExp.cpp.

760{
761 ASSERTL2(i >= 0 && i <= 8, "edge id is out of range");
762
763 if (i == 0 || i == 2)
764 {
765 return GetBasisNumModes(0);
766 }
767 else if (i == 1 || i == 3 || i == 8)
768 {
769 return GetBasisNumModes(1);
770 }
771 else
772 {
773 return GetBasisNumModes(2);
774 }
775}
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:265
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:169

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

Referenced by v_GetEdgeInteriorToElementMap().

◆ v_GetElmtTraceToTraceMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1244 of file StdPrismExp.cpp.

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

References ASSERTL0, ASSERTL1, Nektar::StdRegions::eDir1BwdDir1_Dir2BwdDir2, Nektar::StdRegions::eDir1BwdDir1_Dir2FwdDir2, Nektar::StdRegions::eDir1BwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1BwdDir2_Dir2FwdDir1, Nektar::StdRegions::eDir1FwdDir1_Dir2BwdDir2, Nektar::StdRegions::eDir1FwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1FwdDir2_Dir2FwdDir1, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::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)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1045 of file StdPrismExp.cpp.

1046{
1049 "BasisType is not a boundary interior form");
1052 "BasisType is not a boundary interior form");
1055 "BasisType is not a boundary interior form");
1056
1057 int P = m_base[0]->GetNumModes() - 1, p;
1058 int Q = m_base[1]->GetNumModes() - 1, q;
1059 int R = m_base[2]->GetNumModes() - 1, r;
1060
1061 int nIntCoeffs = m_ncoeffs - NumBndryCoeffs();
1062
1063 if (outarray.size() != nIntCoeffs)
1064 {
1065 outarray = Array<OneD, unsigned int>(nIntCoeffs);
1066 }
1067
1068 int idx = 0;
1069
1070 // Loop over all interior modes.
1071 for (p = 2; p <= P; ++p)
1072 {
1073 for (q = 2; q <= Q; ++q)
1074 {
1075 for (r = 1; r <= R - p; ++r)
1076 {
1077 outarray[idx++] = GetMode(p, q, r);
1078 }
1079 }
1080 }
1081}

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, CellMLToNektar.cellml_metadata::p, and Nektar::UnitTests::q().

◆ v_GetNedges()

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

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 786 of file StdPrismExp.cpp.

787{
788 return 9;
789}

◆ v_GetNtraces()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 791 of file StdPrismExp.cpp.

792{
793 return 5;
794}

◆ v_GetNverts()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 781 of file StdPrismExp.cpp.

782{
783 return 6;
784}

◆ v_GetTraceBasisKey()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 923 of file StdPrismExp.cpp.

925{
926 ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
927 ASSERTL2(k >= 0 && k <= 1, "basis key id is out of range");
928
929 switch (i)
930 {
931 case 0:
932 {
934 m_base[k]->GetNumPoints(),
935 m_base[k]->GetNumModes());
936 }
937 case 2:
938 case 4:
939 {
940 return EvaluateQuadFaceBasisKey(k, m_base[k + 1]->GetBasisType(),
941 m_base[k + 1]->GetNumPoints(),
942 m_base[k + 1]->GetNumModes());
943 }
944 case 1:
945 case 3:
946 {
947 return EvaluateTriFaceBasisKey(k, m_base[2 * k]->GetBasisType(),
948 m_base[2 * k]->GetNumPoints(),
949 m_base[2 * k]->GetNumModes());
950 }
951 break;
952 }
953
954 // Should never get here.
956}
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 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1143 of file StdPrismExp.cpp.

1145{
1147 "Method only implemented if BasisType is identical"
1148 "in x and y directions");
1151 "Method only implemented for Modified_A BasisType"
1152 "(x and y direction) and Modified_B BasisType (z "
1153 "direction)");
1154 int p, q, r, idx = 0;
1155 int P = 0, Q = 0;
1156
1157 switch (fid)
1158 {
1159 case 0:
1160 P = m_base[0]->GetNumModes();
1161 Q = m_base[1]->GetNumModes();
1162 break;
1163 case 1:
1164 case 3:
1165 P = m_base[0]->GetNumModes();
1166 Q = m_base[2]->GetNumModes();
1167 break;
1168 case 2:
1169 case 4:
1170 P = m_base[1]->GetNumModes();
1171 Q = m_base[2]->GetNumModes();
1172 break;
1173 default:
1174 ASSERTL0(false, "fid must be between 0 and 4");
1175 }
1176
1177 if (maparray.size() != P * Q)
1178 {
1179 maparray = Array<OneD, unsigned int>(P * Q);
1180 }
1181
1182 // Set up ordering inside each 2D face. Also for triangular faces,
1183 // populate signarray.
1184 switch (fid)
1185 {
1186 case 0: // Bottom quad
1187 for (q = 0; q < Q; ++q)
1188 {
1189 for (p = 0; p < P; ++p)
1190 {
1191 maparray[q * P + p] = GetMode(p, q, 0);
1192 }
1193 }
1194 break;
1195 case 1: // Left triangle
1196 for (p = 0; p < P; ++p)
1197 {
1198 for (r = 0; r < Q - p; ++r)
1199 {
1200 maparray[idx++] = GetMode(p, 0, r);
1201 }
1202 }
1203 break;
1204 case 2: // Slanted quad
1205 for (q = 0; q < P; ++q)
1206 {
1207 maparray[q] = GetMode(1, q, 0);
1208 }
1209 for (q = 0; q < P; ++q)
1210 {
1211 maparray[P + q] = GetMode(0, q, 1);
1212 }
1213 for (r = 1; r < Q - 1; ++r)
1214 {
1215 for (q = 0; q < P; ++q)
1216 {
1217 maparray[(r + 1) * P + q] = GetMode(1, q, r);
1218 }
1219 }
1220 break;
1221 case 3: // Right triangle
1222 for (p = 0; p < P; ++p)
1223 {
1224 for (r = 0; r < Q - p; ++r)
1225 {
1226 maparray[idx++] = GetMode(p, 1, r);
1227 }
1228 }
1229 break;
1230 case 4: // Rear quad
1231 for (r = 0; r < Q; ++r)
1232 {
1233 for (q = 0; q < P; ++q)
1234 {
1235 maparray[r * P + q] = GetMode(0, q, r);
1236 }
1237 }
1238 break;
1239 default:
1240 ASSERTL0(false, "Face to element map unavailable.");
1241 }
1242}

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

◆ v_GetTraceInteriorToElementMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1632 of file StdPrismExp.cpp.

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

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, Nektar::UnitTests::q(), and v_GetTraceIntNcoeffs().

◆ v_GetTraceIntNcoeffs()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 863 of file StdPrismExp.cpp.

864{
865 ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
866
867 int Pi = GetBasisNumModes(0) - 2;
868 int Qi = GetBasisNumModes(1) - 2;
869 int Ri = GetBasisNumModes(2) - 2;
870
871 if (i == 0)
872 {
873 return Pi * Qi;
874 }
875 else if (i == 1 || i == 3)
876 {
877 return Pi * (2 * Ri - Pi - 1) / 2;
878 }
879 else
880 {
881 return Qi * Ri;
882 }
883}

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

Referenced by v_GetTraceInteriorToElementMap().

◆ v_GetTraceNcoeffs()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 845 of file StdPrismExp.cpp.

846{
847 ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
848 if (i == 0)
849 {
850 return GetBasisNumModes(0) * GetBasisNumModes(1);
851 }
852 else if (i == 1 || i == 3)
853 {
854 int P = GetBasisNumModes(0) - 1, Q = GetBasisNumModes(2) - 1;
855 return Q + 1 + (P * (1 + 2 * Q - P)) / 2;
856 }
857 else
858 {
859 return GetBasisNumModes(1) * GetBasisNumModes(2);
860 }
861}

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 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 721 of file StdPrismExp.cpp.

723{
724 int nummodes[3] = {m_base[0]->GetNumModes(), m_base[1]->GetNumModes(),
725 m_base[2]->GetNumModes()};
726 switch (fid)
727 {
728 // base quad
729 case 0:
730 {
731 numModes0 = nummodes[0];
732 numModes1 = nummodes[1];
733 }
734 break;
735 // front and back quad
736 case 2:
737 case 4:
738 {
739 numModes0 = nummodes[1];
740 numModes1 = nummodes[2];
741 }
742 break;
743 // triangles
744 case 1:
745 case 3:
746 {
747 numModes0 = nummodes[0];
748 numModes1 = nummodes[2];
749 }
750 break;
751 }
752
753 if (faceOrient >= eDir1FwdDir2_Dir2FwdDir1)
754 {
755 std::swap(numModes0, numModes1);
756 }
757}

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

◆ v_GetTraceNumPoints()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 885 of file StdPrismExp.cpp.

886{
887 ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
888
889 if (i == 0)
890 {
891 return m_base[0]->GetNumPoints() * m_base[1]->GetNumPoints();
892 }
893 else if (i == 1 || i == 3)
894 {
895 return m_base[0]->GetNumPoints() * m_base[2]->GetNumPoints();
896 }
897 else
898 {
899 return m_base[1]->GetNumPoints() * m_base[2]->GetNumPoints();
900 }
901}

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

◆ v_GetTracePointsKey()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 903 of file StdPrismExp.cpp.

905{
906 ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
907 ASSERTL2(j == 0 || j == 1, "face direction is out of range");
908
909 if (i == 0)
910 {
911 return m_base[j]->GetPointsKey();
912 }
913 else if (i == 1 || i == 3)
914 {
915 return m_base[2 * j]->GetPointsKey();
916 }
917 else
918 {
919 return m_base[j + 1]->GetPointsKey();
920 }
921}

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

◆ v_GetVertexMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 980 of file StdPrismExp.cpp.

981{
985 "Mapping not defined for this type of basis");
986
987 int l = 0;
988
989 if (useCoeffPacking == true) // follow packing of coefficients i.e q,r,p
990 {
991 switch (vId)
992 {
993 case 0:
994 l = GetMode(0, 0, 0);
995 break;
996 case 1:
997 l = GetMode(0, 0, 1);
998 break;
999 case 2:
1000 l = GetMode(0, 1, 0);
1001 break;
1002 case 3:
1003 l = GetMode(0, 1, 1);
1004 break;
1005 case 4:
1006 l = GetMode(1, 0, 0);
1007 break;
1008 case 5:
1009 l = GetMode(1, 1, 0);
1010 break;
1011 default:
1012 ASSERTL0(false, "local vertex id must be between 0 and 5");
1013 }
1014 }
1015 else
1016 {
1017 switch (vId)
1018 {
1019 case 0:
1020 l = GetMode(0, 0, 0);
1021 break;
1022 case 1:
1023 l = GetMode(1, 0, 0);
1024 break;
1025 case 2:
1026 l = GetMode(1, 1, 0);
1027 break;
1028 case 3:
1029 l = GetMode(0, 1, 0);
1030 break;
1031 case 4:
1032 l = GetMode(0, 0, 1);
1033 break;
1034 case 5:
1035 l = GetMode(0, 1, 1);
1036 break;
1037 default:
1038 ASSERTL0(false, "local vertex id must be between 0 and 5");
1039 }
1040 }
1041
1042 return l;
1043}

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 
)
overrideprotectedvirtual

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

Definition at line 368 of file StdPrismExp.cpp.

370{
373 "Basis[1] is not a general tensor type");
374
377 "Basis[2] is not a general tensor type");
378
379 if (m_base[0]->Collocation() && m_base[1]->Collocation())
380 {
381 MultiplyByQuadratureMetric(inarray, outarray);
382 }
383 else
384 {
385 StdPrismExp::v_IProductWRTBase_SumFac(inarray, outarray);
386 }
387}
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:723
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override

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_SumFac()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PrismExp.

Definition at line 389 of file StdPrismExp.cpp.

392{
393 int nquad1 = m_base[1]->GetNumPoints();
394 int nquad2 = m_base[2]->GetNumPoints();
395 int order0 = m_base[0]->GetNumModes();
396 int order1 = m_base[1]->GetNumModes();
397
398 Array<OneD, NekDouble> wsp(order0 * nquad2 * (nquad1 + order1));
399
400 if (multiplybyweights)
401 {
402 Array<OneD, NekDouble> tmp(inarray.size());
403
404 MultiplyByQuadratureMetric(inarray, tmp);
406 m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
407 tmp, outarray, wsp, true, true, true);
408 }
409 else
410 {
412 m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
413 inarray, outarray, wsp, true, true, true);
414 }
415}
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 
)
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion3D.

Definition at line 417 of file StdPrismExp.cpp.

426{
427 // Interior prism implementation based on Spen's book page
428 // 119. and 608.
429 const int nquad0 = m_base[0]->GetNumPoints();
430 const int nquad1 = m_base[1]->GetNumPoints();
431 const int nquad2 = m_base[2]->GetNumPoints();
432 const int order0 = m_base[0]->GetNumModes();
433 const int order1 = m_base[1]->GetNumModes();
434 const int order2 = m_base[2]->GetNumModes();
435
436 int i, mode;
437
438 ASSERTL1(wsp.size() >= nquad1 * nquad2 * order0 + nquad2 * order0 * order1,
439 "Insufficient workspace size");
440
441 Array<OneD, NekDouble> tmp0 = wsp;
442 Array<OneD, NekDouble> tmp1 = wsp + nquad1 * nquad2 * order0;
443
444 // Inner product with respect to the '0' direction
445 Blas::Dgemm('T', 'N', nquad1 * nquad2, order0, nquad0, 1.0, inarray.get(),
446 nquad0, base0.get(), nquad0, 0.0, tmp0.get(), nquad1 * nquad2);
447
448 // Inner product with respect to the '1' direction
449 Blas::Dgemm('T', 'N', nquad2 * order0, order1, nquad1, 1.0, tmp0.get(),
450 nquad1, base1.get(), nquad1, 0.0, tmp1.get(), nquad2 * order0);
451
452 // Inner product with respect to the '2' direction
453 for (mode = i = 0; i < order0; ++i)
454 {
455 Blas::Dgemm('T', 'N', order2 - i, order1, nquad2, 1.0,
456 base2.get() + mode * nquad2, nquad2,
457 tmp1.get() + i * nquad2, nquad2 * order0, 0.0,
458 outarray.get() + mode * order1, order2 - i);
459 mode += order2 - i;
460 }
461
462 // Fix top singular vertices; performs phi_{0,q,1} +=
463 // phi_1(xi_1)*phi_q(xi_2)*phi_{01}*phi_r(xi_2).
465 {
466 for (i = 0; i < order1; ++i)
467 {
468 mode = GetMode(0, i, 1);
469 outarray[mode] +=
470 Blas::Ddot(nquad2, base2.get() + nquad2, 1,
471 tmp1.get() + i * order0 * nquad2 + nquad2, 1);
472 }
473 }
474}
static double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = .
Definition: Blas.hpp:163

References ASSERTL1, Blas::Ddot(), Blas::Dgemm(), 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 
)
overrideprotectedvirtual

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

Definition at line 480 of file StdPrismExp.cpp.

483{
484 v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
485}
void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override

References v_IProductWRTDerivBase_SumFac().

◆ v_IProductWRTDerivBase_SumFac()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PrismExp.

Definition at line 487 of file StdPrismExp.cpp.

490{
491 ASSERTL0(dir >= 0 && dir <= 2, "input dir is out of range");
492
493 int i;
494 int order0 = m_base[0]->GetNumModes();
495 int order1 = m_base[1]->GetNumModes();
496 int nquad0 = m_base[0]->GetNumPoints();
497 int nquad1 = m_base[1]->GetNumPoints();
498 int nquad2 = m_base[2]->GetNumPoints();
499
500 const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
501 const Array<OneD, const NekDouble> &z2 = m_base[2]->GetZ();
502 Array<OneD, NekDouble> gfac0(nquad0);
503 Array<OneD, NekDouble> gfac2(nquad2);
504 Array<OneD, NekDouble> tmp0(nquad0 * nquad1 * nquad2);
505 Array<OneD, NekDouble> wsp(order0 * nquad2 * (nquad1 + order1));
506
507 // set up geometric factor: (1+z0)/2
508 for (i = 0; i < nquad0; ++i)
509 {
510 gfac0[i] = 0.5 * (1 + z0[i]);
511 }
512
513 // Set up geometric factor: 2/(1-z2)
514 for (i = 0; i < nquad2; ++i)
515 {
516 gfac2[i] = 2.0 / (1 - z2[i]);
517 }
518
519 // Scale first derivative term by gfac2.
520 if (dir != 1)
521 {
522 for (i = 0; i < nquad2; ++i)
523 {
524 Vmath::Smul(nquad0 * nquad1, gfac2[i],
525 &inarray[0] + i * nquad0 * nquad1, 1,
526 &tmp0[0] + i * nquad0 * nquad1, 1);
527 }
528 MultiplyByQuadratureMetric(tmp0, tmp0);
529 }
530
531 switch (dir)
532 {
533 case 0:
534 {
536 m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
537 m_base[2]->GetBdata(), tmp0, outarray, wsp, true, true, true);
538 break;
539 }
540 case 1:
541 {
542 MultiplyByQuadratureMetric(inarray, tmp0);
544 m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
545 m_base[2]->GetBdata(), tmp0, outarray, wsp, true, true, true);
546 break;
547 }
548
549 case 2:
550 {
551 Array<OneD, NekDouble> tmp1(m_ncoeffs);
552
553 // Scale eta_1 derivative with gfac0.
554 for (i = 0; i < nquad1 * nquad2; ++i)
555 {
556 Vmath::Vmul(nquad0, &gfac0[0], 1, &tmp0[0] + i * nquad0, 1,
557 &tmp0[0] + i * nquad0, 1);
558 }
559
561 m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
562 m_base[2]->GetBdata(), tmp0, tmp1, wsp, true, true, true);
563
564 MultiplyByQuadratureMetric(inarray, tmp0);
566 m_base[0]->GetBdata(), m_base[1]->GetBdata(),
567 m_base[2]->GetDbdata(), tmp0, outarray, wsp, true, true, true);
568
569 Vmath::Vadd(m_ncoeffs, &tmp1[0], 1, &outarray[0], 1, &outarray[0],
570 1);
571 break;
572 }
573 }
574}
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.hpp:72
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180

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 ( ) const
overrideprotectedvirtual

◆ v_LocCollapsedToLocCoord()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 600 of file StdPrismExp.cpp.

602{
603 xi[0] = (1.0 + eta[0]) * (1.0 - eta[2]) * 0.5 - 1.0;
604 xi[1] = eta[1];
605 xi[2] = eta[2];
606}

◆ v_LocCoordToLocCollapsed()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 580 of file StdPrismExp.cpp.

582{
583 NekDouble d2 = 1.0 - xi[2];
584 if (fabs(d2) < NekConstants::kNekZeroTol)
585 {
586 if (d2 >= 0.)
587 {
589 }
590 else
591 {
593 }
594 }
595 eta[2] = xi[2]; // eta_z = xi_z
596 eta[1] = xi[1]; // eta_y = xi_y
597 eta[0] = 2.0 * (1.0 + xi[0]) / d2 - 1.0;
598}
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 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1946 of file StdPrismExp.cpp.

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

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

◆ v_NumBndryCoeffs()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 805 of file StdPrismExp.cpp.

806{
809 "BasisType is not a boundary interior form");
812 "BasisType is not a boundary interior form");
815 "BasisType is not a boundary interior form");
816
817 int P = m_base[0]->GetNumModes();
818 int Q = m_base[1]->GetNumModes();
819 int R = m_base[2]->GetNumModes();
820
822}
int getNumberOfBndCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:290

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
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 824 of file StdPrismExp.cpp.

825{
828 "BasisType is not a boundary interior form");
831 "BasisType is not a boundary interior form");
834 "BasisType is not a boundary interior form");
835
836 int P = m_base[0]->GetNumModes() - 1;
837 int Q = m_base[1]->GetNumModes() - 1;
838 int R = m_base[2]->GetNumModes() - 1;
839
840 return (P + 1) * (Q + 1) // 1 rect. face on base
841 + 2 * (Q + 1) * (R + 1) // other 2 rect. faces
842 + 2 * (R + 1) + P * (1 + 2 * R - P); // 2 tri. faces
843}

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 
)
overrideprotectedvirtual

Calculate the derivative of the physical points.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PrismExp.

Definition at line 75 of file StdPrismExp.cpp.

79{
80 int Qx = m_base[0]->GetNumPoints();
81 int Qy = m_base[1]->GetNumPoints();
82 int Qz = m_base[2]->GetNumPoints();
83 int Qtot = Qx * Qy * Qz;
84
85 Array<OneD, NekDouble> dEta_bar1(Qtot, 0.0);
86
87 Array<OneD, const NekDouble> eta_x, eta_z;
88 eta_x = m_base[0]->GetZ();
89 eta_z = m_base[2]->GetZ();
90
91 int i, k;
92
93 bool Do_1 = (out_dxi1.size() > 0) ? true : false;
94 bool Do_3 = (out_dxi3.size() > 0) ? true : false;
95
96 // out_dXi2 is just a tensor derivative so is just passed through
97 if (Do_3)
98 {
99 PhysTensorDeriv(u_physical, dEta_bar1, out_dxi2, out_dxi3);
100 }
101 else if (Do_1)
102 {
103 PhysTensorDeriv(u_physical, dEta_bar1, out_dxi2, NullNekDouble1DArray);
104 }
105 else // case if just require 2nd direction
106 {
107 PhysTensorDeriv(u_physical, NullNekDouble1DArray, out_dxi2,
109 }
110
111 if (Do_1)
112 {
113 for (k = 0; k < Qz; ++k)
114 {
115 Vmath::Smul(Qx * Qy, 2.0 / (1.0 - eta_z[k]),
116 &dEta_bar1[0] + k * Qx * Qy, 1,
117 &out_dxi1[0] + k * Qx * Qy, 1);
118 }
119 }
120
121 if (Do_3)
122 {
123 // divide dEta_Bar1 by (1-eta_z)
124 for (k = 0; k < Qz; ++k)
125 {
126 Vmath::Smul(Qx * Qy, 1.0 / (1.0 - eta_z[k]),
127 &dEta_bar1[0] + k * Qx * Qy, 1,
128 &dEta_bar1[0] + k * Qx * Qy, 1);
129 }
130
131 // Multiply dEta_Bar1 by (1+eta_x) and add ot out_dxi3
132 for (i = 0; i < Qx; ++i)
133 {
134 Vmath::Svtvp(Qz * Qy, 1.0 + eta_x[i], &dEta_bar1[0] + i, Qx,
135 &out_dxi3[0] + i, Qx, &out_dxi3[0] + i, Qx);
136 }
137 }
138}
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.hpp:396

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 
)
overrideprotectedvirtual

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

See also
StdRegions::StdExpansion::PhysDeriv

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 140 of file StdPrismExp.cpp.

143{
144 switch (dir)
145 {
146 case 0:
147 {
148 v_PhysDeriv(inarray, outarray, NullNekDouble1DArray,
150 break;
151 }
152
153 case 1:
154 {
155 v_PhysDeriv(inarray, NullNekDouble1DArray, outarray,
157 break;
158 }
159
160 case 2:
161 {
163 outarray);
164 break;
165 }
166
167 default:
168 {
169 ASSERTL1(false, "input dir is out of range");
170 }
171 break;
172 }
173}
void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2) override
Calculate the derivative of the physical points.
Definition: StdPrismExp.cpp:75

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

◆ v_PhysEvaluate()

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

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Reimplemented in Nektar::LocalRegions::PrismExp.

Definition at line 635 of file StdPrismExp.cpp.

639{
640 // Collapse coordinates
641 Array<OneD, NekDouble> coll(3, 0.0);
642 LocCoordToLocCollapsed(coord, coll);
643
644 // If near singularity do the old interpolation matrix method
645 // @TODO: Dave thinks there might be a way in the Barycentric to
646 // mathematically remove this singularity?
647 if ((1 - coll[2]) < 1e-5)
648 {
649 int totPoints = GetTotPoints();
650 Array<OneD, NekDouble> EphysDeriv0(totPoints), EphysDeriv1(totPoints),
651 EphysDeriv2(totPoints);
652 PhysDeriv(inarray, EphysDeriv0, EphysDeriv1, EphysDeriv2);
653
654 Array<OneD, DNekMatSharedPtr> I(3);
655 I[0] = GetBase()[0]->GetI(coll);
656 I[1] = GetBase()[1]->GetI(coll + 1);
657 I[2] = GetBase()[2]->GetI(coll + 2);
658
659 firstOrderDerivs[0] = PhysEvaluate(I, EphysDeriv0);
660 firstOrderDerivs[1] = PhysEvaluate(I, EphysDeriv1);
661 firstOrderDerivs[2] = PhysEvaluate(I, EphysDeriv2);
662 return PhysEvaluate(I, inarray);
663 }
664
665 NekDouble val = BaryTensorDeriv(coll, inarray, firstOrderDerivs);
666
667 NekDouble dEta_bar1 = firstOrderDerivs[0];
668
669 NekDouble fac = 2.0 / (1.0 - coll[2]);
670 firstOrderDerivs[0] = fac * dEta_bar1;
671
672 // divide dEta_Bar1 by (1-eta_z)
673 fac = 1.0 / (1.0 - coll[2]);
674 dEta_bar1 = fac * dEta_bar1;
675
676 // Multiply dEta_Bar1 by (1+eta_x) and add ot out_dxi3
677 fac = 1.0 + coll[0];
678 firstOrderDerivs[2] += fac * dEta_bar1;
679
680 return val;
681}
NekDouble BaryTensorDeriv(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:134
const Array< OneD, const LibUtilities::BasisSharedPtr > & GetBase() const
This function gets the shared point to basis.
Definition: StdExpansion.h:100
NekDouble PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
This function evaluates the expansion at a single (arbitrary) point of the domain.
Definition: StdExpansion.h:919
void PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Definition: StdExpansion.h:849

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

◆ v_PhysEvaluateBasis()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 690 of file StdPrismExp.cpp.

692{
693 Array<OneD, NekDouble> coll(3);
694 LocCoordToLocCollapsed(coords, coll);
695
696 const int nm1 = m_base[1]->GetNumModes();
697 const int nm2 = m_base[2]->GetNumModes();
698 const int b = 2 * nm2 + 1;
699
700 const int mode0 = floor(0.5 * (b - sqrt(b * b - 8.0 * mode / nm1)));
701 const int tmp =
702 mode - nm1 * (mode0 * (nm2 - 1) + 1 - (mode0 - 2) * (mode0 - 1) / 2);
703 const int mode1 = tmp / (nm2 - mode0);
704 const int mode2 = tmp % (nm2 - mode0);
705
706 if (mode0 == 0 && mode2 == 1 &&
708 {
709 // handle collapsed top edge to remove mode0 terms
710 return StdExpansion::BaryEvaluateBasis<1>(coll[1], mode1) *
711 StdExpansion::BaryEvaluateBasis<2>(coll[2], mode2);
712 }
713 else
714 {
715 return StdExpansion::BaryEvaluateBasis<0>(coll[0], mode0) *
716 StdExpansion::BaryEvaluateBasis<1>(coll[1], mode1) *
717 StdExpansion::BaryEvaluateBasis<2>(coll[2], mode2);
718 }
719}
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

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 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2140 of file StdPrismExp.cpp.

2143{
2144 int nquad0 = m_base[0]->GetNumPoints();
2145 int nquad1 = m_base[1]->GetNumPoints();
2146 int nquad2 = m_base[2]->GetNumPoints();
2147 int nqtot = nquad0 * nquad1 * nquad2;
2148 int nmodes0 = m_base[0]->GetNumModes();
2149 int nmodes1 = m_base[1]->GetNumModes();
2150 int nmodes2 = m_base[2]->GetNumModes();
2151 int numMax = nmodes0;
2152
2153 Array<OneD, NekDouble> coeff(m_ncoeffs);
2154 Array<OneD, NekDouble> coeff_tmp1(m_ncoeffs, 0.0);
2155 Array<OneD, NekDouble> phys_tmp(nqtot, 0.0);
2156 Array<OneD, NekDouble> tmp, tmp2, tmp3, tmp4;
2157
2158 const LibUtilities::PointsKey Pkey0 = m_base[0]->GetPointsKey();
2159 const LibUtilities::PointsKey Pkey1 = m_base[1]->GetPointsKey();
2160 const LibUtilities::PointsKey Pkey2 = m_base[2]->GetPointsKey();
2161
2162 LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A, nmodes0, Pkey0);
2163 LibUtilities::BasisKey bortho1(LibUtilities::eOrtho_A, nmodes1, Pkey1);
2164 LibUtilities::BasisKey bortho2(LibUtilities::eOrtho_B, nmodes2, Pkey2);
2165
2166 int cnt = 0;
2167 int u = 0;
2168 int i = 0;
2170
2172 bortho0, bortho1, bortho2);
2173
2174 BwdTrans(inarray, phys_tmp);
2175 OrthoPrismExp->FwdTrans(phys_tmp, coeff);
2176
2177 // filtering
2178 for (u = 0; u < numMin; ++u)
2179 {
2180 for (i = 0; i < numMin; ++i)
2181 {
2182 Vmath::Vcopy(numMin - u, tmp = coeff + cnt, 1,
2183 tmp2 = coeff_tmp1 + cnt, 1);
2184 cnt += numMax - u;
2185 }
2186
2187 for (i = numMin; i < numMax; ++i)
2188 {
2189 cnt += numMax - u;
2190 }
2191 }
2192
2193 OrthoPrismExp->BwdTrans(coeff_tmp1, phys_tmp);
2194 StdPrismExp::FwdTrans(phys_tmp, outarray);
2195}
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Backward transformation from coefficient space to physical space.
Definition: StdExpansion.h:424
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Forward transformation from physical space to coefficient space.
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:42
std::shared_ptr< StdPrismExp > StdPrismExpSharedPtr
Definition: StdPrismExp.h:218

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 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 175 of file StdPrismExp.cpp.

179{
180 StdPrismExp::v_PhysDeriv(inarray, out_d0, out_d1, out_d2);
181}

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 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 183 of file StdPrismExp.cpp.

186{
187 StdPrismExp::v_PhysDeriv(dir, inarray, outarray);
188}

References v_PhysDeriv().

◆ v_SVVLaplacianFilter()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PrismExp.

Definition at line 2002 of file StdPrismExp.cpp.

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

References Nektar::StdRegions::StdExpansion::BwdTrans(), Nektar::StdRegions::StdMatrixKey::ConstFactorExists(), Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::eFactorSVVDGKerDiffCoeff, Nektar::StdRegions::eFactorSVVDiffCoeff, Nektar::StdRegions::eFactorSVVPowerKerDiffCoeff, Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::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.