Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Protected Member Functions | Private Attributes | List of all members
Nektar::StdRegions::StdPyrExp Class Reference

#include <StdPyrExp.h>

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

Public Member Functions

 StdPyrExp ()
 StdPyrExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
 StdPyrExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc, NekDouble *coeffs, NekDouble *phys)
 StdPyrExp (const StdPyrExp &T)
 ~StdPyrExp ()
int GetTetMode (int I, int J, int K)
 Number tetrahedral modes in r-direction. Much the same as StdTetExp::GetTetMode but slightly simplified since we know that the polynomial order is the same in each direction.
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion3D
 StdExpansion3D ()
 StdExpansion3D (int numcoeffs, const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
 StdExpansion3D (const StdExpansion3D &T)
virtual ~StdExpansion3D ()
void PhysTensorDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d1, Array< OneD, NekDouble > &outarray_d2, Array< OneD, NekDouble > &outarray_d3)
 Calculate the 3D derivative in the local tensor/collapsed coordinate at the physical points.
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)
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion
 StdExpansion ()
 Default Constructor.
 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.
 StdExpansion (const StdExpansion &T)
 Copy Constructor.
virtual ~StdExpansion ()
 Destructor.
int GetNumBases () const
 This function returns the number of 1D bases used in the expansion.
const Array< OneD, const
LibUtilities::BasisSharedPtr > & 
GetBase () const
 This function gets the shared point to basis.
const
LibUtilities::BasisSharedPtr
GetBasis (int dir) const
 This function gets the shared point to basis in the dir direction.
int GetNcoeffs (void) const
 This function returns the total number of coefficients used in the expansion.
int GetTotPoints () const
 This function returns the total number of quadrature points used in the element.
LibUtilities::BasisType GetBasisType (const int dir) const
 This function returns the type of basis used in the dir direction.
int GetBasisNumModes (const int dir) const
 This function returns the number of expansion modes in the dir direction.
int EvalBasisNumModesMax (void) const
 This function returns the maximum number of expansion modes over all local directions.
LibUtilities::PointsType GetPointsType (const int dir) const
 This function returns the type of quadrature points used in the dir direction.
int GetNumPoints (const int dir) const
 This function returns the number of quadrature points in the dir direction.
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.
int GetNverts () const
 This function returns the number of vertices of the expansion domain.
int GetNedges () const
 This function returns the number of edges of the expansion domain.
int GetEdgeNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th edge.
int GetTotalEdgeIntNcoeffs () const
int GetEdgeNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th edge.
int DetCartesianDirOfEdge (const int edge)
const LibUtilities::BasisKey DetEdgeBasisKey (const int i) const
const LibUtilities::BasisKey DetFaceBasisKey (const int i, const int k) const
int GetFaceNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th face.
int GetFaceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th face.
int GetFaceIntNcoeffs (const int i) const
int GetTotalFaceIntNcoeffs () const
LibUtilities::PointsKey GetFacePointsKey (const int i, const int j) const
int NumBndryCoeffs (void) const
int NumDGBndryCoeffs (void) const
LibUtilities::BasisType GetEdgeBasisType (const int i) const
 This function returns the type of expansion basis on the i-th edge.
int GetNfaces () const
 This function returns the number of faces of the expansion domain.
LibUtilities::ShapeType DetShapeType () const
 This function returns the shape of the expansion domain.
int GetShapeDimension () const
bool IsBoundaryInteriorExpansion ()
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space.
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space.
void FwdTrans_BndConstrained (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.
void FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 This function fills the array outarray with the mode-th mode of the expansion.
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
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)
int GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id.
void SetElmtId (const int id)
 Set the element id of this expansion when used in a list by returning value of m_elmt_id.
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
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
DNekMatSharedPtr GetStdMatrix (const StdMatrixKey &mkey)
DNekBlkMatSharedPtr GetStdStaticCondMatrix (const StdMatrixKey &mkey)
IndexMapValuesSharedPtr GetIndexMap (const IndexMapKey &ikey)
const Array< OneD, const
NekDouble > & 
GetPhysNormals (void)
void SetPhysNormals (Array< OneD, const NekDouble > &normal)
virtual void SetUpPhysNormals (const int edge)
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const 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)
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
void DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
StdRegions::Orientation GetFaceOrient (int face)
StdRegions::Orientation GetEorient (int edge)
StdRegions::Orientation GetPorient (int point)
StdRegions::Orientation GetCartesianEorient (int edge)
void SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
void SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
int CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
void ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs)
NekDouble StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
void AddEdgeNormBoundaryInt (const int edge, boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
void AddEdgeNormBoundaryInt (const int edge, boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
void AddEdgeNormBoundaryBiInt (const int edge, boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
void AddFaceNormBoundaryInt (const int face, boost::shared_ptr< StdExpansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
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 GetEdgeInteriorMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
void GetFaceInteriorMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
void GetEdgeToElementMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
void GetFaceToElementMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int nummodesA=-1, int nummodesB=-1)
void GetEdgePhysVals (const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Extract the physical values along edge edge from inarray into outarray following the local edge orientation and point distribution defined by defined in EdgeExp.
void GetEdgePhysVals (const int edge, const boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void GetTracePhysVals (const int edge, const boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void GetVertexPhysVals (const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray)
void GetEdgeInterpVals (const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void GetEdgeQFactors (const int edge, Array< OneD, NekDouble > &outarray)
 Extract the metric factors to compute the contravariant fluxes along edge edge and stores them into outarray following the local edge orientation (i.e. anticlockwise convention).
void GetFacePhysVals (const int face, const boost::shared_ptr< StdExpansion > &FaceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=eNoOrientation)
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}$
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 LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void WeakDirectionalDerivMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void MassLevelCurvatureMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LinearAdvectionDiffusionReactionMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
void HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
DNekMatSharedPtr GenMatrix (const StdMatrixKey &mkey)
void PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void PhysDeriv_s (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds)
void PhysDeriv_n (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn)
void PhysDirectionalDeriv (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &outarray)
void StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
void StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void AddRobinMassMatrix (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
void AddRobinEdgeContribution (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, Array< OneD, NekDouble > &coeffs)
void DGDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, boost::shared_ptr< StdExpansion > > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &coeffs, 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.
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.
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta.
const boost::shared_ptr
< SpatialDomains::GeomFactors > & 
GetMetricInfo (void) const
virtual int v_GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id.
virtual const Array< OneD,
const NekDouble > & 
v_GetPhysNormals (void)
virtual void v_SetPhysNormals (Array< OneD, const NekDouble > &normal)
virtual void v_SetUpPhysNormals (const int edge)
virtual void v_ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmode_offset, NekDouble *coeffs)
 Unpack data from input file assuming it comes from the same expansion type.
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 DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
virtual void v_DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
virtual StdRegions::Orientation v_GetFaceOrient (int face)
virtual StdRegions::Orientation v_GetEorient (int edge)
virtual StdRegions::Orientation v_GetCartesianEorient (int edge)
virtual StdRegions::Orientation v_GetPorient (int point)
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.
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.
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.
const NormalVectorGetEdgeNormal (const int edge) const
void ComputeEdgeNormal (const int edge)
void NegateEdgeNormal (const int edge)
bool EdgeNormalNegated (const int edge)
void ComputeFaceNormal (const int face)
void NegateFaceNormal (const int face)
void ComputeVertexNormal (const int vertex)
const NormalVectorGetFaceNormal (const int face) const
const NormalVectorGetVertexNormal (const int vertex) const
const NormalVectorGetSurfaceNormal (const int id) const
const LibUtilities::PointsKeyVector GetPointsKeys () const
Array< OneD, unsigned int > GetEdgeInverseBoundaryMap (int eid)
Array< OneD, unsigned int > GetFaceInverseBoundaryMap (int fid, StdRegions::Orientation faceOrient=eNoOrientation)
DNekMatSharedPtr BuildInverseTransformationMatrix (const DNekScalMatSharedPtr &m_transformationmatrix)
template<class T >
boost::shared_ptr< T > as ()

Protected Member Functions

virtual void v_PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
 Calculate the derivative of the physical points.
virtual void v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Calculate the derivative of the physical points in a given direction.
virtual void v_StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
virtual void v_StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Backward transformation is evaluated at the quadrature points.
virtual void v_BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
virtual void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Forward transform from physical quadrature space stored in inarray and evaluate the expansion coefficients and store in outarray.
virtual void v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Inner product of inarray over region with respect to the expansion basis m_base[0]->GetBdata(),m_base[1]->GetBdata(), m_base[2]->GetBdata() and return in outarray.
virtual void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
virtual void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
virtual void v_GetCoords (Array< OneD, NekDouble > &xi_x, Array< OneD, NekDouble > &xi_y, Array< OneD, NekDouble > &xi_z)
virtual void v_FillMode (const int mode, Array< OneD, NekDouble > &outarray)
virtual int v_GetNverts () const
virtual int v_GetNedges () const
virtual int v_GetNfaces () const
virtual LibUtilities::ShapeType v_DetShapeType () const
virtual int v_NumBndryCoeffs () const
virtual int v_GetEdgeNcoeffs (const int i) const
virtual int v_GetFaceNcoeffs (const int i) const
virtual int v_GetFaceIntNcoeffs (const int i) const
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
virtual LibUtilities::BasisType v_GetEdgeBasisType (const int i) const
virtual void v_GetFaceToElementMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int nummodesA=-1, int nummodesB=-1)
virtual int v_GetVertexMap (int localVertexId, bool useCoeffPacking=false)
virtual void v_GetEdgeInteriorMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
virtual void v_GetFaceInteriorMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
virtual void v_GetInteriorMap (Array< OneD, unsigned int > &outarray)
virtual void v_GetBoundaryMap (Array< OneD, unsigned int > &outarray)
virtual DNekMatSharedPtr v_GenMatrix (const StdMatrixKey &mkey)
virtual DNekMatSharedPtr v_CreateStdMatrix (const StdMatrixKey &mkey)
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion3D
virtual NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain.
virtual NekDouble v_PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
virtual void v_LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
virtual void v_HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
virtual NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray)
 Integrates the specified function over the domain.
virtual void v_NegateFaceNormal (const int face)
- 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.
IndexMapValuesSharedPtr CreateIndexMap (const IndexMapKey &ikey)
 Create an IndexMap which contains mapping information linking any specific element shape with either its boundaries, edges, faces, verteces, etc.
void BwdTrans_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void IProductWRTBase_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 GeneralMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void MassMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
void LaplacianMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LaplacianMatrixOp_MatFree (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void WeakDerivMatrixOp_MatFree (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void WeakDirectionalDerivMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void MassLevelCurvatureMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LinearAdvectionDiffusionReactionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
void HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void HelmholtzMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
virtual NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
virtual void v_AddEdgeNormBoundaryInt (const int edge, boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
virtual void v_AddEdgeNormBoundaryInt (const int edge, boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
virtual void v_AddFaceNormBoundaryInt (const int face, boost::shared_ptr< StdExpansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)

Private Attributes

map< Mode, unsigned int, cmpopm_map
map< int, map< int, map< int,
pair< int, int > > > > 
m_idxMap

Additional Inherited Members

- Protected Attributes inherited from Nektar::StdRegions::StdExpansion3D
std::map< int, NormalVectorm_faceNormals
std::map< int, bool > m_negatedNormals

Detailed Description

Definition at line 82 of file StdPyrExp.h.

Constructor & Destructor Documentation

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

Definition at line 44 of file StdPyrExp.cpp.

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

Definition at line 48 of file StdPyrExp.cpp.

References ASSERTL0, Nektar::StdRegions::find(), Nektar::LibUtilities::BasisKey::GetNumModes(), GetTetMode(), Nektar::iterator, m_idxMap, m_map, and Nektar::StdRegions::StdExpansion::m_ncoeffs.

Ba.GetNumModes(),
Bb.GetNumModes(),
Bc.GetNumModes()),
3, Ba, Bb, Bc),
Ba.GetNumModes(),
Bb.GetNumModes(),
Bc.GetNumModes()),
Ba, Bb, Bc)
{
if (Ba.GetNumModes() > Bc.GetNumModes())
{
ASSERTL0(false, "order in 'a' direction is higher "
"than order in 'c' direction");
}
if (Bb.GetNumModes() > Bc.GetNumModes())
{
ASSERTL0(false, "order in 'b' direction is higher "
"than order in 'c' direction");
}
// Set up mode mapping which takes 0\leq i\leq N_coeffs -> (p,q,r)
// of the 3D tensor product
const int P = Ba.GetNumModes() - 1;
const int Q = Bb.GetNumModes() - 1;
const int R = Bc.GetNumModes() - 1;
int cnt = 0;
// Vertices
m_map[Mode(0, 0, 0, 0)] = cnt++;
m_map[Mode(1, 0, 0, 0)] = cnt++;
m_map[Mode(1, 1, 0, 0)] = cnt++;
m_map[Mode(0, 1, 0, 0)] = cnt++;
m_map[Mode(0, 0, 1, 1)] = cnt++;
// Edge 0
for (int i = 2; i <= P; ++i)
{
m_map[Mode(i, 0, 0, GetTetMode(i, 0, 0))] = cnt++;
}
// Edge 1
for (int i = 2; i <= Q; ++i)
{
m_map[Mode(1, i, 0, GetTetMode(0, i, 0))] = cnt++;
}
// Edge 2
for (int i = 2; i <= P; ++i)
{
m_map[Mode(i, 1, 0, GetTetMode(i, 0, 0))] = cnt++;
}
// Edge 3
for (int i = 2; i <= Q; ++i)
{
m_map[Mode(0, i, 0, GetTetMode(0, i, 0))] = cnt++;
}
// Edge 4
for (int i = 2; i <= R; ++i)
{
m_map[Mode(0, 0, i, i)] = cnt++;
}
// Edge 5
for (int i = 2; i <= R; ++i)
{
m_map[Mode(1, 0, i, i)] = cnt++;
}
// Edge 6
for (int i = 2; i <= R; ++i)
{
m_map[Mode(1, 1, i, i)] = cnt++;
}
// Edge 7
for (int i = 2; i <= R; ++i)
{
m_map[Mode(0, 1, i, i)] = cnt++;
}
// Face 0 - TODO check this
for (int j = 2; j <= Q; ++j)
{
for (int i = 2; i <= P; ++i)
{
m_map[Mode(i, j, 0, GetTetMode((i-2+j-2) % (Q-1) + 2, 0, 0))] = cnt++;
}
}
// Face 1
for (int i = 2; i <= P; ++i)
{
for (int j = 1; j <= R-i; ++j)
{
m_map[Mode(i, 0, j, GetTetMode(i, 0, j))] = cnt++;
}
}
// Face 2
for (int i = 2; i <= Q; ++i)
{
for (int j = 1; j <= R-i; ++j)
{
m_map[Mode(1, i, j, GetTetMode(0, i, j))] = cnt++;
}
}
// Face 3
for (int i = 2; i <= P; ++i)
{
for (int j = 1; j <= R-i; ++j)
{
m_map[Mode(i, 1, j, GetTetMode(i, 0, j))] = cnt++;
}
}
// Face 4
for (int i = 2; i <= Q; ++i)
{
for (int j = 1; j <= R-i; ++j)
{
m_map[Mode(0, i, j, GetTetMode(0, i, j))] = cnt++;
}
}
// Interior (tetrahedral modes)
for (int i = 2; i <= P+1; ++i)
{
for (int j = 1; j <= Q-i+1; ++j)
{
for (int k = 1; k <= R-i-j+1; ++k)
{
// need to go to j+1-th mode in the 'b' direction to
// select correct modified_a mode
m_map[Mode(i, j+1, k, GetTetMode(i-1, j, k))] = cnt++;
}
}
}
"Duplicate coefficient entries in map");
for (it = m_map.begin(); it != m_map.end(); ++it)
{
const int p = it->first.get<0>();
const int q = it->first.get<1>();
const int r = it->first.get<2>();
const int rp = it->first.get<3>();
if (m_idxMap.find(p) == m_idxMap.end())
{
m_idxMap[p] = map<int, map<int, pair<int, int> > >();
}
if (m_idxMap[p].find(q) == m_idxMap[p].end())
{
m_idxMap[p][q] = map<int, pair<int, int> >();
}
if (m_idxMap[p][q].find(r) == m_idxMap[p][q].end())
{
m_idxMap[p][q][r] = pair<int, int>(it->second, rp);
}
}
}
Nektar::StdRegions::StdPyrExp::StdPyrExp ( const LibUtilities::BasisKey Ba,
const LibUtilities::BasisKey Bb,
const LibUtilities::BasisKey Bc,
NekDouble coeffs,
NekDouble phys 
)
Nektar::StdRegions::StdPyrExp::StdPyrExp ( const StdPyrExp T)

Definition at line 221 of file StdPyrExp.cpp.

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

Definition at line 229 of file StdPyrExp.cpp.

{
}

Member Function Documentation

int Nektar::StdRegions::StdPyrExp::GetTetMode ( int  I,
int  J,
int  K 
)

Number tetrahedral modes in r-direction. Much the same as StdTetExp::GetTetMode but slightly simplified since we know that the polynomial order is the same in each direction.

Definition at line 1789 of file StdPyrExp.cpp.

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

Referenced by StdPyrExp().

{
const int R = m_base[2]->GetNumModes();
int i, j, cnt = 0;
for (i = 0; i < I; ++i)
{
cnt += (R-i)*(R-i+1)/2;
}
i = R-I;
for (j = 0; j < J; ++j)
{
cnt += i;
i--;
}
return cnt + K;
}
void Nektar::StdRegions::StdPyrExp::v_BwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Backward transformation is evaluated at the quadrature points.

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

Backward transformation is three dimensional tensorial expansion

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

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 367 of file StdPyrExp.cpp.

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

Referenced by v_FillMode().

{
if (m_base[0]->Collocation() &&
m_base[1]->Collocation() &&
m_base[2]->Collocation())
{
inarray, 1, outarray, 1);
}
else
{
StdPyrExp::v_BwdTrans_SumFac(inarray,outarray);
}
}
void Nektar::StdRegions::StdPyrExp::v_BwdTrans_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Sum-factorisation implementation of the BwdTrans operation.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 388 of file StdPyrExp.cpp.

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

Referenced by v_BwdTrans().

{
Array<OneD, NekDouble> wsp;
m_base[1]->GetBdata(),
m_base[2]->GetBdata(),
inarray,outarray,wsp,
true,true,true);
}
void Nektar::StdRegions::StdPyrExp::v_BwdTrans_SumFacKernel ( const Array< OneD, const NekDouble > &  base0,
const Array< OneD, const NekDouble > &  base1,
const Array< OneD, const NekDouble > &  base2,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
Array< OneD, NekDouble > &  wsp,
bool  doCheckCollDir0,
bool  doCheckCollDir1,
bool  doCheckCollDir2 
)
protectedvirtual

Implements Nektar::StdRegions::StdExpansion3D.

Definition at line 400 of file StdPyrExp.cpp.

References Nektar::iterator, Nektar::StdRegions::StdExpansion::m_base, and m_idxMap.

Referenced by v_BwdTrans_SumFac().

{
const int Qx = m_base[0]->GetNumPoints();
const int Qy = m_base[1]->GetNumPoints();
const int Qz = m_base[2]->GetNumPoints();
const NekDouble *bx = base0.get();
const NekDouble *by = base1.get();
const NekDouble *bz = base2.get();
// Need to count coeffs for storage...
map<int, map<int, map<int, pair<int, int> > > >::iterator it_p;
map<int, map<int, pair<int, int> > > ::iterator it_q;
map<int, pair<int, int> > ::iterator it_r;
int pqcnt = 0;
for (it_p = m_idxMap.begin(); it_p != m_idxMap.end(); ++it_p)
{
for (it_q = it_p->second.begin(); it_q != it_p->second.end(); ++it_q)
{
pqcnt++;
}
}
Array<OneD, NekDouble> fpq(pqcnt);
Array<OneD, NekDouble> fp (m_base[0]->GetNumModes());
int i ,j, k, s = 0, cnt = 0, cnt2 = 0;
for (k = 0; k < Qz; ++k)
{
NekDouble bz1 = bz[k+Qz];
cnt = 0;
for (it_p = m_idxMap.begin(); it_p != m_idxMap.end(); ++it_p)
{
for (it_q = it_p->second.begin(); it_q != it_p->second.end(); ++it_q)
{
NekDouble sum = 0.0;
for (it_r = it_q->second.begin(); it_r != it_q->second.end(); ++it_r)
{
sum += inarray[it_r->second.first] * bz[k + Qz*it_r->second.second];
}
fpq[cnt++] = sum;
}
}
for (j = 0; j < Qy; ++j)
{
NekDouble by0 = bz1*by[j];
NekDouble by1 = bz1*by[j+Qy];
cnt = cnt2 = 0;
for (it_p = m_idxMap.begin(); it_p != m_idxMap.end(); ++it_p)
{
NekDouble sum = 0.0;
for (it_q = it_p->second.begin(); it_q != it_p->second.end(); ++it_q)
{
sum += by[j + Qy*it_q->first] * fpq[cnt++];
}
fp[cnt2++] = sum;
}
for (i = 0; i < Qx; ++i, ++s)
{
cnt2 = 0;
NekDouble sum = 0.0;
for (it_p = m_idxMap.begin(); it_p != m_idxMap.end(); ++it_p)
{
sum += bx[i + Qx*it_p->first] * fp[cnt2++];
}
sum += inarray[4]*(by1*(bx[i] + bx[i+Qx]) + by0*bx[i+Qx]);
outarray[s] = sum;
}
}
}
}
int Nektar::StdRegions::StdPyrExp::v_CalcNumberOfCoefficients ( const std::vector< unsigned int > &  nummodes,
int &  modes_offset 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 961 of file StdPyrExp.cpp.

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

{
nummodes[modes_offset],
nummodes[modes_offset+1],
nummodes[modes_offset+2]);
modes_offset += 3;
return nmodes;
}
DNekMatSharedPtr Nektar::StdRegions::StdPyrExp::v_CreateStdMatrix ( const StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PyrExp.

Definition at line 1779 of file StdPyrExp.cpp.

References v_GenMatrix().

{
return v_GenMatrix(mkey);
}
LibUtilities::ShapeType Nektar::StdRegions::StdPyrExp::v_DetShapeType ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 876 of file StdPyrExp.cpp.

References Nektar::LibUtilities::ePyramid.

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 849 of file StdPyrExp.cpp.

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

{
Array<OneD, NekDouble> tmp(m_ncoeffs, 0.0);
tmp[mode] = 1.0;
v_BwdTrans(tmp, outarray);
}
void Nektar::StdRegions::StdPyrExp::v_FwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

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

Inputs:

  • inarray: array of physical quadrature points to be transformed

Outputs:

  • outarray: updated array of expansion coefficients.

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PyrExp.

Definition at line 499 of file StdPyrExp.cpp.

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

{
v_IProductWRTBase(inarray,outarray);
// get Mass matrix inverse
StdMatrixKey masskey(eInvMass,DetShapeType(),*this);
DNekMatSharedPtr matsys = GetStdMatrix(masskey);
// copy inarray in case inarray == outarray
DNekVec in (m_ncoeffs, outarray);
DNekVec out(m_ncoeffs, outarray, eWrapper);
out = (*matsys)*in;
}
DNekMatSharedPtr Nektar::StdRegions::StdPyrExp::v_GenMatrix ( const StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PyrExp.

Definition at line 1774 of file StdPyrExp.cpp.

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

Referenced by v_CreateStdMatrix().

{
return CreateGeneralMatrix(mkey);
}
void Nektar::StdRegions::StdPyrExp::v_GetBoundaryMap ( Array< OneD, unsigned int > &  outarray)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1750 of file StdPyrExp.cpp.

References ASSERTL1, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_C, Nektar::StdRegions::StdExpansion::GetBasisType(), and v_NumBndryCoeffs().

{
"BasisType is not a boundary interior form");
"BasisType is not a boundary interior form");
"BasisType is not a boundary interior form");
int idx = 0, nBndry = v_NumBndryCoeffs();
for (idx = 0; idx < nBndry; ++idx)
{
maparray[idx] = idx;
}
}
void Nektar::StdRegions::StdPyrExp::v_GetCoords ( Array< OneD, NekDouble > &  xi_x,
Array< OneD, NekDouble > &  xi_y,
Array< OneD, NekDouble > &  xi_z 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PyrExp.

Definition at line 821 of file StdPyrExp.cpp.

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

{
Array<OneD, const NekDouble> etaBar_x = m_base[0]->GetZ();
Array<OneD, const NekDouble> eta_y = m_base[1]->GetZ();
Array<OneD, const NekDouble> eta_z = m_base[2]->GetZ();
int Qx = GetNumPoints(0);
int Qy = GetNumPoints(1);
int Qz = GetNumPoints(2);
// Convert collapsed coordinates into cartesian coordinates: eta --> xi
for (int k = 0; k < Qz; ++k )
{
for (int j = 0; j < Qy; ++j)
{
for (int i = 0; i < Qx; ++i)
{
int s = i + Qx*(j + Qy*k);
xi_z[s] = eta_z[k];
xi_y[s] = (1.0 + eta_y[j]) * (1.0 - eta_z[k]) / 2.0 - 1.0;
xi_x[s] = (1.0 + etaBar_x[i]) * (1.0 - eta_z[k]) / 2.0 - 1.0;
}
}
}
}
LibUtilities::BasisType Nektar::StdRegions::StdPyrExp::v_GetEdgeBasisType ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 974 of file StdPyrExp.cpp.

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

{
ASSERTL2(i >= 0 && i <= 7, "edge id is out of range");
if (i == 0 || i == 2)
{
return GetBasisType(0);
}
else if (i == 1 || i == 3)
{
return GetBasisType(1);
}
else
{
return GetBasisType(2);
}
}
void Nektar::StdRegions::StdPyrExp::v_GetEdgeInteriorMap ( const int  eid,
const Orientation  edgeOrient,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1480 of file StdPyrExp.cpp.

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

{
int i;
bool signChange;
const int P = m_base[0]->GetNumModes() - 2;
const int Q = m_base[1]->GetNumModes() - 2;
const int R = m_base[2]->GetNumModes() - 2;
const int nEdgeIntCoeffs = v_GetEdgeNcoeffs(eid) - 2;
if (maparray.num_elements() != nEdgeIntCoeffs)
{
maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
}
if(signarray.num_elements() != nEdgeIntCoeffs)
{
signarray = Array<OneD, int>(nEdgeIntCoeffs,1);
}
else
{
fill(signarray.get(), signarray.get()+nEdgeIntCoeffs, 1);
}
// If edge is oriented backwards, change sign of modes which have
// degree 2n+1, n >= 1.
signChange = edgeOrient == eBackwards;
int offset = 5;
switch (eid)
{
case 0:
break;
case 1:
offset += P;
break;
case 2:
offset += P+Q;
break;
case 3:
offset += 2*P+Q;
break;
case 4:
offset += 2*(P+Q);
break;
case 5:
offset += 2*(P+Q)+R;
break;
case 6:
offset += 2*(P+Q+R);
break;
case 7:
offset += 2*(P+Q)+3*R;
break;
default:
ASSERTL0(false, "Edge not defined.");
break;
}
for (i = 0; i < nEdgeIntCoeffs; ++i)
{
maparray[i] = offset + i;
}
if (signChange)
{
for (i = 1; i < nEdgeIntCoeffs; i += 2)
{
signarray[i] = -1;
}
}
}
int Nektar::StdRegions::StdPyrExp::v_GetEdgeNcoeffs ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 901 of file StdPyrExp.cpp.

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

Referenced by v_GetEdgeInteriorMap().

{
ASSERTL2(i >= 0 && i <= 7, "edge id is out of range");
if (i == 0 || i == 2)
{
return GetBasisNumModes(0);
}
else if (i == 1 || i == 3)
{
return GetBasisNumModes(1);
}
else
{
return GetBasisNumModes(2);
}
}
void Nektar::StdRegions::StdPyrExp::v_GetFaceInteriorMap ( const int  fid,
const Orientation  faceOrient,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1557 of file StdPyrExp.cpp.

References ASSERTL0, Nektar::StdRegions::StdExpansion::m_base, and v_GetFaceIntNcoeffs().

{
const int P = m_base[0]->GetNumModes() - 1;
const int Q = m_base[1]->GetNumModes() - 1;
const int R = m_base[2]->GetNumModes() - 1;
const int nFaceIntCoeffs = v_GetFaceIntNcoeffs(fid);
int p, q, r, idx = 0;
int nummodesA;
int nummodesB;
int i, j;
if (maparray.num_elements() != nFaceIntCoeffs)
{
maparray = Array<OneD, unsigned int>(nFaceIntCoeffs);
}
if (signarray.num_elements() != nFaceIntCoeffs)
{
signarray = Array<OneD, int>(nFaceIntCoeffs, 1);
}
else
{
fill(signarray.get(), signarray.get()+nFaceIntCoeffs, 1);
}
// Set up an array indexing for quad faces, since the ordering may
// need to be transposed depending on orientation.
Array<OneD, int> arrayindx(nFaceIntCoeffs);
if (fid == 0)
{
nummodesA = P-1;
nummodesB = Q-1;
for (i = 0; i < nummodesB; i++)
{
for (j = 0; j < nummodesA; j++)
{
if (faceOrient < 9)
{
arrayindx[i*nummodesA+j] = i*nummodesA+j;
}
else
{
arrayindx[i*nummodesA+j] = j*nummodesB+i;
}
}
}
}
int offset = 5 + 2*(P-1) + 2*(Q-1) + 4*(R-1);
for (i = 0; i < fid; ++i)
{
offset += v_GetFaceIntNcoeffs(i);
}
switch (fid)
{
case 0:
for (q = 2; q <= Q; ++q)
{
for (p = 2; p <= P; ++p)
{
maparray[arrayindx[(q-2)*nummodesA+(p-2)]]
= offset + (q-2)*nummodesA+(p-2);
}
}
break;
case 1:
case 3:
for (p = 2; p <= P; ++p)
{
for (r = 1; r <= R-p; ++r, ++idx)
{
if ((int)faceOrient == 7)
{
signarray[idx] = p % 2 ? -1 : 1;
}
maparray[idx] = offset + idx;
}
}
break;
case 2:
case 4:
for (q = 2; q <= Q; ++q)
{
for (r = 1; r <= R-q; ++r, ++idx)
{
if ((int)faceOrient == 7)
{
signarray[idx] = q % 2 ? -1 : 1;
}
maparray[idx] = offset + idx;
}
}
break;
default:
ASSERTL0(false, "Face interior map not available.");
}
// Triangular faces are processed in the above switch loop; for
// remaining quad faces, set up orientation if necessary.
if (fid > 0)
{
return;
}
if (faceOrient == 6 || faceOrient == 8 ||
faceOrient == 11 || faceOrient == 12)
{
if (faceOrient < 9)
{
for (i = 1; i < nummodesB; i += 2)
{
for (j = 0; j < nummodesA; j++)
{
signarray[arrayindx[i*nummodesA+j]] *= -1;
}
}
}
else
{
for (i = 0; i < nummodesB; i++)
{
for (j = 1; j < nummodesA; j += 2)
{
signarray[arrayindx[i*nummodesA+j]] *= -1;
}
}
}
}
if (faceOrient == 7 || faceOrient == 8 ||
faceOrient == 10 || faceOrient == 12)
{
if (faceOrient < 9)
{
for (i = 0; i < nummodesB; i++)
{
for (j = 1; j < nummodesA; j += 2)
{
signarray[arrayindx[i*nummodesA+j]] *= -1;
}
}
}
else
{
for (i = 1; i < nummodesB; i += 2)
{
for (j = 0; j < nummodesA; j++)
{
signarray[arrayindx[i*nummodesA+j]] *= -1;
}
}
}
}
}
int Nektar::StdRegions::StdPyrExp::v_GetFaceIntNcoeffs ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 939 of file StdPyrExp.cpp.

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

Referenced by v_GetFaceInteriorMap(), and v_GetFaceToElementMap().

{
ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
int P = m_base[0]->GetNumModes()-1;
int Q = m_base[1]->GetNumModes()-1;
int R = m_base[2]->GetNumModes()-1;
if (i == 0)
{
return (P-1)*(Q-1);
}
else if (i == 1 || i == 3)
{
return (P-1) * (2*(R-1) - (P-1) - 1) / 2;
}
else
{
return (Q-1) * (2*(R-1) - (Q-1) - 1) / 2;
}
}
int Nektar::StdRegions::StdPyrExp::v_GetFaceNcoeffs ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 919 of file StdPyrExp.cpp.

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

{
ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
if (i == 0)
{
}
else if (i == 1 || i == 3)
{
int P = GetBasisNumModes(0)-1, Q = GetBasisNumModes(2)-1;
return Q+1 + (P*(1 + 2*Q - P))/2;
}
else
{
int P = GetBasisNumModes(1)-1, Q = GetBasisNumModes(2)-1;
return Q+1 + (P*(1 + 2*Q - P))/2;
}
}
void Nektar::StdRegions::StdPyrExp::v_GetFaceToElementMap ( const int  fid,
const Orientation  faceOrient,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray,
int  nummodesA = -1,
int  nummodesB = -1 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 996 of file StdPyrExp.cpp.

References ASSERTL0, ASSERTL1, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_C, Nektar::StdRegions::StdExpansion::GetEdgeBasisType(), Nektar::StdRegions::StdExpansion::GetFaceNcoeffs(), Nektar::StdRegions::StdExpansion::m_base, and v_GetFaceIntNcoeffs().

{
"Method only implemented if BasisType is identical"
"in x and y directions");
"Method only implemented for Modified_A BasisType"
"(x and y direction) and Modified_C BasisType (z "
"direction)");
int i, j, p, q, r, nFaceCoeffs;
int order0 = m_base[0]->GetNumModes();
int order1 = m_base[1]->GetNumModes();
int order2 = m_base[2]->GetNumModes();
if (nummodesA == -1)
{
switch (fid)
{
case 0:
nummodesA = m_base[0]->GetNumModes();
nummodesB = m_base[1]->GetNumModes();
break;
case 1:
case 3:
nummodesA = m_base[0]->GetNumModes();
nummodesB = m_base[2]->GetNumModes();
break;
case 2:
case 4:
nummodesA = m_base[1]->GetNumModes();
nummodesB = m_base[2]->GetNumModes();
break;
}
nFaceCoeffs = GetFaceNcoeffs(fid);
}
else if (fid > 0)
{
nFaceCoeffs = nummodesB + (nummodesA-1)*(1+2*(nummodesB-1)-(nummodesA-1))/2;
}
else
{
nFaceCoeffs = nummodesA*nummodesB;
}
// Allocate the map array and sign array; set sign array to ones (+)
if (maparray.num_elements() != nFaceCoeffs)
{
maparray = Array<OneD, unsigned int>(nFaceCoeffs);
}
if (signarray.num_elements() != nFaceCoeffs)
{
signarray = Array<OneD, int>(nFaceCoeffs,1);
}
else
{
fill(signarray.get(), signarray.get() + nFaceCoeffs, 1);
}
// Set up an array indexing for quads, since the ordering may need
// to be transposed.
Array<OneD, int> arrayindx(nFaceCoeffs,-1);
if (fid == 0)
{
for (i = 0; i < nummodesB; i++)
{
for (j = 0; j < nummodesA; j++)
{
if (faceOrient < 9)
{
arrayindx[i*nummodesA+j] = i*nummodesA+j;
}
else
{
arrayindx[i*nummodesA+j] = j*nummodesB+i;
}
}
}
}
// Set up ordering inside each 2D face. Also for triangular faces,
// populate signarray.
int cnt = 0, cnt2;
switch (fid)
{
case 0: // Bottom quad
// Fill in vertices
maparray[arrayindx[0]] = 0;
maparray[arrayindx[1]] = 1;
maparray[arrayindx[nummodesA+1]] = 2;
maparray[arrayindx[nummodesA]] = 3;
// Edge 0
cnt = 5;
for (p = 2; p < nummodesA; ++p)
{
maparray[arrayindx[p]] = p-2 + cnt;
}
// Edge 1
cnt += nummodesA-2;
for (q = 2; q < nummodesB; ++q)
{
maparray[arrayindx[q*nummodesA+1]] = q-2 + cnt;
}
// Edge 2
cnt += nummodesB-2;
for (p = 2; p < nummodesA; ++p)
{
maparray[arrayindx[nummodesA+p]] = p-2 + cnt;
}
// Edge 3
cnt += nummodesA-2;
for (q = 2; q < nummodesB; ++q)
{
maparray[arrayindx[q*nummodesA]] = q-2 + cnt;
}
// Interior
cnt += nummodesB-2 + 4*(nummodesA-2);
for (q = 2; q < nummodesB; ++q)
{
for (p = 2; p < nummodesA; ++p)
{
maparray[arrayindx[q*nummodesA+p]] = cnt + (q-2)*nummodesA+(p-2);
}
}
break;
case 1: // Left triangle
// Vertices
maparray[0] = 0;
maparray[1] = 4;
maparray[nummodesB] = 1;
// Edge 0 (pyramid edge 0)
cnt = 5;
q = 2*nummodesB-1;
for (p = 2; p < nummodesA; q += nummodesB-p, ++p)
{
maparray[q] = cnt++;
if ((int)faceOrient == 7)
{
signarray[q] = p % 2 ? -1 : 1;
}
}
// Edge 1 (pyramid edge 5)
cnt = 5 + 2*(order0-2) + 2*(order1-2) + (order2-2);
for (q = 2; q < nummodesB; ++q)
{
maparray[q] = cnt++;
if ((int)faceOrient == 7)
{
signarray[q] = q % 2 ? -1 : 1;
}
}
// Edge 2 (pyramid edge 4)
cnt = 5 + 2*(order0-2) + 2*(order1-2);
for (q = 2; q < nummodesB; ++q)
{
maparray[nummodesB+q-1] = cnt++;
if ((int)faceOrient == 7)
{
signarray[nummodesB+q-1] = q % 2 ? -1 : 1;
}
}
// Interior
cnt = 5 + 2*(order0-2) + 2*(order1-2) + 4*(order2-2)
cnt2 = 2*nummodesB + 1;
for (p = 2; p < nummodesA; ++p)
{
for (r = 2; r < nummodesB-p; ++r)
{
maparray[cnt2] = cnt++;
if ((int)faceOrient == 7 && p > 1)
{
signarray[cnt2++] = p % 2 ? -1 : 1;
}
}
cnt2++;
}
break;
case 2:
// Vertices
maparray[0] = 1;
maparray[1] = 4;
maparray[nummodesB] = 2;
// Edge 0 (pyramid edge 1)
cnt = 5 + (order0-2);
q = 2*nummodesB-1;
for (p = 2; p < nummodesA; q += nummodesB-p, ++p)
{
maparray[q] = cnt++;
if ((int)faceOrient == 7)
{
signarray[q] = p % 2 ? -1 : 1;
}
}
// Edge 1 (pyramid edge 6)
cnt = 5 + 2*(order0-2) + 2*(order1-2) + 2*(order2-2);
for (q = 2; q < nummodesB; ++q)
{
maparray[q] = cnt++;
if ((int)faceOrient == 7)
{
signarray[q] = q % 2 ? -1 : 1;
}
}
// Edge 2 (pyramid edge 5)
cnt = 5 + 2*(order0-2) + 2*(order1-2) + (order2-2);
for (q = 2; q < nummodesB; ++q)
{
maparray[nummodesB+q-1] = cnt++;
if ((int)faceOrient == 7)
{
signarray[nummodesB+q-1] = q % 2 ? -1 : 1;
}
}
// Interior
cnt = 5 + 2*(order0-2) + 2*(order1-2) + 4*(order2-2)
cnt2 = 2*nummodesB + 1;
for (p = 2; p < nummodesA; ++p)
{
for (r = 2; r < nummodesB-p; ++r)
{
maparray[cnt2] = cnt++;
if ((int)faceOrient == 7 && p > 1)
{
signarray[cnt2++] = p % 2 ? -1 : 1;
}
}
cnt2++;
}
break;
case 3: // Right triangle
// Vertices
maparray[0] = 3;
maparray[1] = 4;
maparray[nummodesB] = 2;
// Edge 0 (pyramid edge 2)
cnt = 5 + (order0-2) + (order1-2);
q = 2*nummodesB-1;
for (p = 2; p < nummodesA; q += nummodesB-p, ++p)
{
maparray[q] = cnt++;
if ((int)faceOrient == 7)
{
signarray[q] = p % 2 ? -1 : 1;
}
}
// Edge 1 (pyramid edge 6)
cnt = 5 + 2*(order0-2) + 2*(order1-2) + 2*(order2-2);
for (q = 2; q < nummodesB; ++q)
{
maparray[q] = cnt++;
if ((int)faceOrient == 7)
{
signarray[q] = q % 2 ? -1 : 1;
}
}
// Edge 2 (pyramid edge 7)
cnt = 5 + 2*(order0-2) + 2*(order1-2) + 3*(order2-2);
for (q = 2; q < nummodesB; ++q)
{
maparray[nummodesB+q-1] = cnt++;
if ((int)faceOrient == 7)
{
signarray[nummodesB+q-1] = q % 2 ? -1 : 1;
}
}
// Interior
cnt = 5 + 2*(order0-2) + 2*(order1-2) + 4*(order2-2)
cnt2 = 2*nummodesB + 1;
for (p = 2; p < nummodesA; ++p)
{
for (r = 2; r < nummodesB-p; ++r)
{
maparray[cnt2] = cnt++;
if ((int)faceOrient == 7 && p > 1)
{
signarray[cnt2++] = p % 2 ? -1 : 1;
}
}
cnt2++;
}
break;
case 4: // Rear quad
// Vertices
maparray[0] = 0;
maparray[1] = 4;
maparray[nummodesB] = 3;
// Edge 0 (pyramid edge 3)
cnt = 5 + 2*(order0-2) + (order1-2);
q = 2*nummodesB-1;
for (p = 2; p < nummodesA; q += nummodesB-p, ++p)
{
maparray[q] = cnt++;
if ((int)faceOrient == 7)
{
signarray[q] = p % 2 ? -1 : 1;
}
}
// Edge 1 (pyramid edge 7)
cnt = 5 + 2*(order0-2) + 2*(order1-2) + 3*(order2-2);
for (q = 2; q < nummodesB; ++q)
{
maparray[q] = cnt++;
if ((int)faceOrient == 7)
{
signarray[q] = q % 2 ? -1 : 1;
}
}
// Edge 2 (pyramid edge 4)
cnt = 5 + 2*(order0-2) + 2*(order1-2);
for (q = 2; q < nummodesB; ++q)
{
maparray[nummodesB+q-1] = cnt++;
if ((int)faceOrient == 7)
{
signarray[nummodesB+q-1] = q % 2 ? -1 : 1;
}
}
// Interior
cnt = 5 + 2*(order0-2) + 2*(order1-2) + 4*(order2-2)
cnt2 = 2*nummodesB + 1;
for (p = 2; p < nummodesA; ++p)
{
for (r = 2; r < nummodesB-p; ++r)
{
maparray[cnt2] = cnt++;
if ((int)faceOrient == 7 && p > 1)
{
signarray[cnt2++] = p % 2 ? -1 : 1;
}
}
cnt2++;
}
break;
default:
ASSERTL0(false, "Face to element map unavailable.");
}
if (fid > 0)
{
// Triangles only have one possible orientation (base
// direction reversed); swap edge modes.
if ((int)faceOrient == 7)
{
swap(maparray[0], maparray[nummodesA]);
for (i = 1; i < nummodesA-1; ++i)
{
swap(maparray[i+1], maparray[nummodesA+i]);
}
}
}
else
{
// The code below is exactly the same as that taken from
// StdHexExp and reverses the 'b' and 'a' directions as
// appropriate (1st and 2nd if statements respectively) in
// quadrilateral faces.
if (faceOrient == 6 || faceOrient == 8 ||
faceOrient == 11 || faceOrient == 12)
{
if (faceOrient < 9)
{
for (i = 3; i < nummodesB; i += 2)
{
for (j = 0; j < nummodesA; j++)
{
signarray[arrayindx[i*nummodesA+j]] *= -1;
}
}
for (i = 0; i < nummodesA; i++)
{
swap(maparray [i], maparray [i+nummodesA]);
swap(signarray[i], signarray[i+nummodesA]);
}
}
else
{
for (i = 0; i < nummodesB; i++)
{
for (j = 3; j < nummodesA; j += 2)
{
signarray[arrayindx[i*nummodesA+j]] *= -1;
}
}
for (i = 0; i < nummodesB; i++)
{
swap (maparray [i], maparray [i+nummodesB]);
swap (signarray[i], signarray[i+nummodesB]);
}
}
}
if (faceOrient == 7 || faceOrient == 8 ||
faceOrient == 10 || faceOrient == 12)
{
if (faceOrient < 9)
{
for (i = 0; i < nummodesB; i++)
{
for (j = 3; j < nummodesA; j += 2)
{
signarray[arrayindx[i*nummodesA+j]] *= -1;
}
}
for(i = 0; i < nummodesB; i++)
{
swap(maparray [i*nummodesA], maparray [i*nummodesA+1]);
swap(signarray[i*nummodesA], signarray[i*nummodesA+1]);
}
}
else
{
for (i = 3; i < nummodesB; i += 2)
{
for (j = 0; j < nummodesA; j++)
{
signarray[arrayindx[i*nummodesA+j]] *= -1;
}
}
for (i = 0; i < nummodesA; i++)
{
swap(maparray [i*nummodesB], maparray [i*nummodesB+1]);
swap(signarray[i*nummodesB], signarray[i*nummodesB+1]);
}
}
}
}
}
void Nektar::StdRegions::StdPyrExp::v_GetInteriorMap ( Array< OneD, unsigned int > &  outarray)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1722 of file StdPyrExp.cpp.

References ASSERTL1, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_C, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::StdRegions::StdExpansion::NumBndryCoeffs(), and v_NumBndryCoeffs().

{
"BasisType is not a boundary interior form");
"BasisType is not a boundary interior form");
"BasisType is not a boundary interior form");
const int nBndCoeffs = v_NumBndryCoeffs();
const int nIntCoeffs = m_ncoeffs - NumBndryCoeffs();
if (outarray.num_elements() != nIntCoeffs)
{
outarray = Array<OneD, unsigned int>(nIntCoeffs);
}
// Loop over all interior modes.
int p, idx = 0;
for (p = nBndCoeffs; p < m_ncoeffs; ++p)
{
outarray[idx++] = p;
}
}
int Nektar::StdRegions::StdPyrExp::v_GetNedges ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 866 of file StdPyrExp.cpp.

{
return 8;
}
int Nektar::StdRegions::StdPyrExp::v_GetNfaces ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 871 of file StdPyrExp.cpp.

{
return 5;
}
int Nektar::StdRegions::StdPyrExp::v_GetNverts ( ) const
protectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 861 of file StdPyrExp.cpp.

{
return 5;
}
int Nektar::StdRegions::StdPyrExp::v_GetVertexMap ( int  localVertexId,
bool  useCoeffPacking = false 
)
protectedvirtual
void Nektar::StdRegions::StdPyrExp::v_IProductWRTBase ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

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

Wrapper call to StdPyrExp::IProductWRTBase

Input:

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

Output:

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

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PyrExp.

Definition at line 534 of file StdPyrExp.cpp.

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

Referenced by v_FwdTrans().

{
if (m_base[0]->Collocation() &&
m_base[1]->Collocation() &&
m_base[2]->Collocation())
{
MultiplyByStdQuadratureMetric(inarray, outarray);
}
else
{
}
}
void Nektar::StdRegions::StdPyrExp::v_IProductWRTBase_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 550 of file StdPyrExp.cpp.

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

Referenced by v_IProductWRTBase().

{
Array<OneD, NekDouble> tmp(inarray.num_elements());
Array<OneD, NekDouble> wsp;
m_base[1]->GetBdata(),
m_base[2]->GetBdata(),
tmp,outarray,wsp,
true,true,true);
}
void Nektar::StdRegions::StdPyrExp::v_IProductWRTBase_SumFacKernel ( const Array< OneD, const NekDouble > &  base0,
const Array< OneD, const NekDouble > &  base1,
const Array< OneD, const NekDouble > &  base2,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
Array< OneD, NekDouble > &  wsp,
bool  doCheckCollDir0,
bool  doCheckCollDir1,
bool  doCheckCollDir2 
)
protectedvirtual

Implements Nektar::StdRegions::StdExpansion3D.

Definition at line 566 of file StdPyrExp.cpp.

References Nektar::iterator, Nektar::StdRegions::StdExpansion::m_base, and m_idxMap.

Referenced by v_IProductWRTBase_SumFac().

{
int i, j, k, s;
int Qx = m_base[0]->GetNumPoints();
int Qy = m_base[1]->GetNumPoints();
int Qz = m_base[2]->GetNumPoints();
const NekDouble *bx = base0.get();
const NekDouble *by = base1.get();
const NekDouble *bz = base2.get();
map<int, map<int, map<int, pair<int, int> > > >::iterator it_p;
map<int, map<int, pair<int, int> > > ::iterator it_q;
map<int, pair<int, int> > ::iterator it_r;
Array<OneD, NekDouble> f (Qy*Qz);
Array<OneD, NekDouble> fb(Qz);
for (it_p = m_idxMap.begin(); it_p != m_idxMap.end(); ++it_p)
{
const int p = it_p->first;
s = 0;
for (k = 0; k < Qz; ++k)
{
for (j = 0; j < Qy; ++j)
{
NekDouble sum = 0.0;
for (i = 0; i < Qx; ++i, ++s)
{
sum += bx[i + Qx*p]*inarray[s];
}
f[j+Qy*k] = sum;
}
}
for (it_q = it_p->second.begin(); it_q != it_p->second.end(); ++it_q)
{
const int q = it_q->first;
for (k = 0; k < Qz; ++k)
{
NekDouble sum = 0.0;
for (j = 0; j < Qy; ++j)
{
sum += by[j + Qy*q]*f[j+Qy*k];
}
fb[k] = sum;
}
for (it_r = it_q->second.begin(); it_r != it_q->second.end(); ++it_r)
{
const int rpqr = it_r->second.second;
NekDouble sum = 0.0;
for (k = 0; k < Qz; ++k)
{
sum += bz[k + Qz*rpqr]*fb[k];
}
outarray[it_r->second.first] = sum;
}
}
}
// Correct for top mode
s = 0;
for (k = 0; k < Qz; ++k)
{
for (j = 0; j < Qy; ++j)
{
for (i = 0; i < Qx; ++i, ++s)
{
outarray[4] += inarray[s] * bz[k+Qz]*(
bx[i+Qx]*by[j+Qy] +
bx[i+Qx]*by[j ] +
bx[i ]*by[j+Qy]);
}
}
}
}
void Nektar::StdRegions::StdPyrExp::v_IProductWRTDerivBase ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 655 of file StdPyrExp.cpp.

References v_IProductWRTDerivBase_SumFac().

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 669 of file StdPyrExp.cpp.

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

Referenced by v_IProductWRTDerivBase().

{
int i;
int nquad0 = m_base[0]->GetNumPoints();
int nquad1 = m_base[1]->GetNumPoints();
int nquad2 = m_base[2]->GetNumPoints();
int nqtot = nquad0*nquad1*nquad2;
Array<OneD, NekDouble> gfac0(nquad0);
Array<OneD, NekDouble> gfac1(nquad1);
Array<OneD, NekDouble> gfac2(nquad2);
Array<OneD, NekDouble> tmp0 (nqtot);
Array<OneD, NekDouble> wsp;
const Array<OneD, const NekDouble>& z0 = m_base[0]->GetZ();
const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
const Array<OneD, const NekDouble>& z2 = m_base[2]->GetZ();
// set up geometric factor: (1+z0)/2
for(i = 0; i < nquad0; ++i)
{
gfac0[i] = 0.5*(1+z0[i]);
}
// set up geometric factor: (1+z1)/2
for(i = 0; i < nquad1; ++i)
{
gfac1[i] = 0.5*(1+z1[i]);
}
// Set up geometric factor: 2/(1-z2)
for(i = 0; i < nquad2; ++i)
{
gfac2[i] = 2.0/(1-z2[i]);
}
// Derivative in first/second direction is always scaled as follows
const int nq01 = nquad0*nquad1;
for(i = 0; i < nquad2; ++i)
{
Vmath::Smul(nq01, gfac2[i],
&inarray[0] + i*nq01, 1,
&tmp0 [0] + i*nq01, 1);
}
switch(dir)
{
case 0:
{
m_base[1]->GetBdata (),
m_base[2]->GetBdata (),
tmp0, outarray, wsp,
false, true, true);
break;
}
case 1:
{
m_base[1]->GetDbdata(),
m_base[2]->GetBdata (),
tmp0, outarray, wsp,
true, false, true);
break;
}
case 2:
{
Array<OneD, NekDouble> tmp3(m_ncoeffs);
Array<OneD, NekDouble> tmp4(m_ncoeffs);
// Scale eta_1 derivative by gfac0
for(i = 0; i < nquad1*nquad2; ++i)
{
Vmath::Vmul(nquad0, tmp0 .get() + i*nquad0, 1,
gfac0.get(), 1,
tmp0 .get() + i*nquad0, 1);
}
m_base[1]->GetBdata(),
m_base[2]->GetBdata(),
tmp0, tmp3, wsp,
false, true, true);
// Scale eta_2 derivative by gfac1*gfac2
for(i = 0; i < nquad2; ++i)
{
Vmath::Smul(nq01, gfac2[i],
&inarray[0] + i*nq01, 1,
&tmp0 [0] + i*nq01, 1);
}
for(i = 0; i < nquad1*nquad2; ++i)
{
Vmath::Smul(nquad0, gfac1[i%nquad1],
&tmp0[0] + i*nquad0, 1,
&tmp0[0] + i*nquad0, 1);
}
m_base[1]->GetDbdata(),
m_base[2]->GetBdata(),
tmp0, tmp4, wsp,
true, false, true);
m_base[1]->GetBdata(),
m_base[2]->GetDbdata(),
tmp0,outarray,wsp,
true, true, false);
Vmath::Vadd(m_ncoeffs,&tmp3[0],1,&outarray[0],1,&outarray[0],1);
Vmath::Vadd(m_ncoeffs,&tmp4[0],1,&outarray[0],1,&outarray[0],1);
break;
}
default:
{
ASSERTL1(false, "input dir is out of range");
break;
}
}
}
void Nektar::StdRegions::StdPyrExp::v_LocCoordToLocCollapsed ( const Array< OneD, const NekDouble > &  xi,
Array< OneD, NekDouble > &  eta 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 801 of file StdPyrExp.cpp.

References Nektar::NekConstants::kNekZeroTol.

{
if (fabs(xi[2]-1.0) < NekConstants::kNekZeroTol)
{
// Very top point of the pyramid
eta[0] = -1.0;
eta[1] = -1.0;
eta[2] = xi[2];
}
else
{
// Below the line-singularity -- Common case
eta[2] = xi[2]; // eta_z = xi_z
eta[1] = 2.0*(1.0 + xi[1])/(1.0 - xi[2]) - 1.0;
eta[0] = 2.0*(1.0 + xi[0])/(1.0 - xi[2]) - 1.0;
}
}
void Nektar::StdRegions::StdPyrExp::v_MultiplyByStdQuadratureMetric ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1808 of file StdPyrExp.cpp.

References Nektar::LibUtilities::eGaussRadauMAlpha2Beta0, Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::StdRegions::StdExpansion::m_base, and Vmath::Vmul().

Referenced by v_IProductWRTBase_SumFac().

{
int i, j;
int nquad0 = m_base[0]->GetNumPoints();
int nquad1 = m_base[1]->GetNumPoints();
int nquad2 = m_base[2]->GetNumPoints();
const Array<OneD, const NekDouble>& w0 = m_base[0]->GetW();
const Array<OneD, const NekDouble>& w1 = m_base[1]->GetW();
const Array<OneD, const NekDouble>& w2 = m_base[2]->GetW();
const Array<OneD, const NekDouble>& z2 = m_base[2]->GetZ();
// Multiply by integration constants in x-direction
for(i = 0; i < nquad1*nquad2; ++i)
{
Vmath::Vmul(nquad0, inarray.get()+i*nquad0, 1,
w0.get(), 1, outarray.get()+i*nquad0,1);
}
// Multiply by integration constants in y-direction
for(j = 0; j < nquad2; ++j)
{
for(i = 0; i < nquad1; ++i)
{
Blas::Dscal(nquad0,w1[i], &outarray[0]+i*nquad0 +
j*nquad0*nquad1,1);
}
}
// Multiply by integration constants in z-direction; need to
// incorporate factor [(1-eta_3)/2]^2 into weights, but only if
// using GLL quadrature points.
switch(m_base[2]->GetPointsType())
{
// (2,0) Jacobi inner product.
for(i = 0; i < nquad2; ++i)
{
Blas::Dscal(nquad0*nquad1, 0.25*w2[i],
&outarray[0]+i*nquad0*nquad1, 1);
}
break;
default:
for(i = 0; i < nquad2; ++i)
{
Blas::Dscal(nquad0*nquad1,0.125*(1-z2[i])*(1-z2[i])*w2[i],
&outarray[0]+i*nquad0*nquad1,1);
}
break;
}
}
int Nektar::StdRegions::StdPyrExp::v_NumBndryCoeffs ( ) const
protectedvirtual
void Nektar::StdRegions::StdPyrExp::v_PhysDeriv ( const Array< OneD, const NekDouble > &  u_physical,
Array< OneD, NekDouble > &  out_dxi1,
Array< OneD, NekDouble > &  out_dxi2,
Array< OneD, NekDouble > &  out_dxi3 
)
protectedvirtual

Calculate the derivative of the physical points.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PyrExp.

Definition at line 251 of file StdPyrExp.cpp.

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

Referenced by v_PhysDeriv(), and v_StdPhysDeriv().

{
// PhysDerivative implementation based on Spen's book page 152.
int Qx = m_base[0]->GetNumPoints();
int Qy = m_base[1]->GetNumPoints();
int Qz = m_base[2]->GetNumPoints();
Array<OneD, NekDouble> dEta_bar1(Qx*Qy*Qz,0.0);
Array<OneD, NekDouble> dXi2 (Qx*Qy*Qz,0.0);
Array<OneD, NekDouble> dEta3 (Qx*Qy*Qz,0.0);
PhysTensorDeriv(u_physical, dEta_bar1, dXi2, dEta3);
Array<OneD, const NekDouble> eta_x, eta_y, eta_z;
eta_x = m_base[0]->GetZ();
eta_y = m_base[1]->GetZ();
eta_z = m_base[2]->GetZ();
int i, j, k, n;
for (k = 0, n = 0; k < Qz; ++k)
{
for (j = 0; j < Qy; ++j)
{
for (i = 0; i < Qx; ++i, ++n)
{
if (out_dxi1.num_elements() > 0)
out_dxi1[n] = 2.0/(1.0 - eta_z[k]) * dEta_bar1[n];
if (out_dxi2.num_elements() > 0)
out_dxi2[n] = 2.0/(1.0 - eta_z[k]) * dXi2[n];
if (out_dxi3.num_elements() > 0)
out_dxi3[n] = (1.0+eta_x[i])/(1.0-eta_z[k])*dEta_bar1[n] +
(1.0+eta_y[j])/(1.0-eta_z[k])*dXi2[n] + dEta3[n];
}
}
}
}
void Nektar::StdRegions::StdPyrExp::v_PhysDeriv ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d0 
)
protectedvirtual

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

See Also
StdRegions::StdExpansion::PhysDeriv

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 292 of file StdPyrExp.cpp.

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

{
switch(dir)
{
case 0:
{
v_PhysDeriv(inarray, outarray, NullNekDouble1DArray,
break;
}
case 1:
{
v_PhysDeriv(inarray, NullNekDouble1DArray, outarray,
break;
}
case 2:
{
break;
}
default:
{
ASSERTL1(false,"input dir is out of range");
}
break;
}
}
void Nektar::StdRegions::StdPyrExp::v_StdPhysDeriv ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d0,
Array< OneD, NekDouble > &  out_d1,
Array< OneD, NekDouble > &  out_d2 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 327 of file StdPyrExp.cpp.

References v_PhysDeriv().

{
StdPyrExp::v_PhysDeriv(inarray, out_d0, out_d1, out_d2);
}
void Nektar::StdRegions::StdPyrExp::v_StdPhysDeriv ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 335 of file StdPyrExp.cpp.

References v_PhysDeriv().

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

Member Data Documentation

map<int, map<int, map<int, pair<int, int> > > > Nektar::StdRegions::StdPyrExp::m_idxMap
private

Definition at line 252 of file StdPyrExp.h.

Referenced by StdPyrExp(), v_BwdTrans_SumFacKernel(), and v_IProductWRTBase_SumFacKernel().

map<Mode, unsigned int, cmpop> Nektar::StdRegions::StdPyrExp::m_map
private

Definition at line 251 of file StdPyrExp.h.

Referenced by StdPyrExp().