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. More...
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion3D
 StdExpansion3D ()
 
 StdExpansion3D (int numcoeffs, const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
 
 StdExpansion3D (const StdExpansion3D &T)
 
virtual ~StdExpansion3D ()
 
void PhysTensorDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d1, Array< OneD, NekDouble > &outarray_d2, Array< OneD, NekDouble > &outarray_d3)
 Calculate the 3D derivative in the local tensor/collapsed coordinate at the physical points. More...
 
void BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
 
void IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
 
- 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::BasisSharedPtr
GetBasis (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 GetNedges () const
 This function returns the number of edges of the expansion domain. More...
 
int GetEdgeNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th edge. More...
 
int GetTotalEdgeIntNcoeffs () const
 
int GetEdgeNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th edge. More...
 
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. More...
 
int GetFaceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th face. More...
 
int GetFaceIntNcoeffs (const int i) const
 
int GetTotalFaceIntNcoeffs () const
 
int GetTraceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th edge/face. More...
 
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. More...
 
const LibUtilities::PointsKey GetNodalPointsKey () const
 This function returns the type of expansion Nodal point type if defined. More...
 
int GetNfaces () const
 This function returns the number of faces of the expansion domain. More...
 
int GetNtrace () 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...
 
boost::shared_ptr< StdExpansionGetStdExp (void) const
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion ()
 
bool IsNodalNonTensorialExp ()
 
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space. More...
 
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space. More...
 
void 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. 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)
 
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)
 
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, 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)
 
StdRegions::Orientation GetForient (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)
 
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, int P=-1)
 
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. More...
 
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). More...
 
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 GetEdgePhysMap (const int edge, Array< OneD, int > &outarray)
 
void GetFacePhysMap (const int face, Array< OneD, int > &outarray)
 
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 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)
 
NekDouble PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta. More...
 
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. More...
 
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. More...
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual void v_DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual StdRegions::Orientation v_GetForient (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. 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 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)
 
bool FaceNormalNegated (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)
 
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 >
boost::shared_ptr< T > as ()
 
void IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 

Protected Member Functions

virtual void v_PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
 Calculate the derivative of the physical points. More...
 
virtual void v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Calculate the derivative of the physical points in a given direction. More...
 
virtual void v_StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
 
virtual void v_StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
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. More...
 
virtual void v_BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
 
virtual void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Forward transform from physical quadrature space stored in inarray and evaluate the expansion coefficients and store in outarray. More...
 
virtual void v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Inner product of inarray over region with respect to the expansion basis m_base[0]->GetBdata(),m_base[1]->GetBdata(), m_base[2]->GetBdata() and return in outarray. More...
 
virtual void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
virtual void v_IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
 
virtual void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
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_GetFaceNumPoints (const int i) const
 
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
virtual const
LibUtilities::BasisKey 
v_DetFaceBasisKey (const int i, const int k) const
 
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. More...
 
virtual NekDouble v_PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 
virtual void v_LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
 
virtual void v_HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
 
virtual NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray)
 Integrates the specified function over the domain. More...
 
virtual void v_NegateFaceNormal (const int face)
 
virtual bool v_FaceNormalNegated (const int face)
 
virtual int v_GetTraceNcoeffs (const int i) const
 
- 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...
 
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. More...
 
void BwdTrans_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void 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)
 

Private Attributes

std::map< Mode, unsigned int,
cmpop
m_map
 
std::map< int, std::map< int,
std::map< int, std::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
 
- Protected Attributes inherited from Nektar::StdRegions::StdExpansion
Array< OneD,
LibUtilities::BasisSharedPtr
m_base
 
int m_elmt_id
 
int m_ncoeffs
 
LibUtilities::NekManager
< StdMatrixKey, DNekMat,
StdMatrixKey::opLess
m_stdMatrixManager
 
LibUtilities::NekManager
< StdMatrixKey, DNekBlkMat,
StdMatrixKey::opLess
m_stdStaticCondMatrixManager
 
LibUtilities::NekManager
< IndexMapKey, IndexMapValues,
IndexMapKey::opLess
m_IndexMapManager
 

Detailed Description

Definition at line 82 of file StdPyrExp.h.

Constructor & Destructor Documentation

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

Definition at line 46 of file StdPyrExp.cpp.

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

Definition at line 50 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.

54  Ba.GetNumModes(),
55  Bb.GetNumModes(),
56  Bc.GetNumModes()),
57  3, Ba, Bb, Bc),
59  Ba.GetNumModes(),
60  Bb.GetNumModes(),
61  Bc.GetNumModes()),
62  Ba, Bb, Bc)
63  {
64  if (Ba.GetNumModes() > Bc.GetNumModes())
65  {
66  ASSERTL0(false, "order in 'a' direction is higher "
67  "than order in 'c' direction");
68  }
69  if (Bb.GetNumModes() > Bc.GetNumModes())
70  {
71  ASSERTL0(false, "order in 'b' direction is higher "
72  "than order in 'c' direction");
73  }
74 
75  // Set up mode mapping which takes 0\leq i\leq N_coeffs -> (p,q,r)
76  // of the 3D tensor product
77  const int P = Ba.GetNumModes() - 1;
78  const int Q = Bb.GetNumModes() - 1;
79  const int R = Bc.GetNumModes() - 1;
80  int cnt = 0;
81 
82  // Vertices
83  m_map[Mode(0, 0, 0, 0)] = cnt++;
84  m_map[Mode(1, 0, 0, 0)] = cnt++;
85  m_map[Mode(1, 1, 0, 0)] = cnt++;
86  m_map[Mode(0, 1, 0, 0)] = cnt++;
87  m_map[Mode(0, 0, 1, 1)] = cnt++;
88 
89  // Edge 0
90  for (int i = 2; i <= P; ++i)
91  {
92  m_map[Mode(i, 0, 0, GetTetMode(i, 0, 0))] = cnt++;
93  }
94 
95  // Edge 1
96  for (int i = 2; i <= Q; ++i)
97  {
98  m_map[Mode(1, i, 0, GetTetMode(0, i, 0))] = cnt++;
99  }
100 
101  // Edge 2
102  for (int i = 2; i <= P; ++i)
103  {
104  m_map[Mode(i, 1, 0, GetTetMode(i, 0, 0))] = cnt++;
105  }
106 
107  // Edge 3
108  for (int i = 2; i <= Q; ++i)
109  {
110  m_map[Mode(0, i, 0, GetTetMode(0, i, 0))] = cnt++;
111  }
112 
113  // Edge 4
114  for (int i = 2; i <= R; ++i)
115  {
116  m_map[Mode(0, 0, i, i)] = cnt++;
117  }
118 
119  // Edge 5
120  for (int i = 2; i <= R; ++i)
121  {
122  m_map[Mode(1, 0, i, i)] = cnt++;
123  }
124 
125  // Edge 6
126  for (int i = 2; i <= R; ++i)
127  {
128  m_map[Mode(1, 1, i, i)] = cnt++;
129  }
130 
131  // Edge 7
132  for (int i = 2; i <= R; ++i)
133  {
134  m_map[Mode(0, 1, i, i)] = cnt++;
135  }
136 
137  // Face 0 - TODO check this
138  for (int j = 2; j <= Q; ++j)
139  {
140  for (int i = 2; i <= P; ++i)
141  {
142  m_map[Mode(i, j, 0, GetTetMode((i-2+j-2) % (Q-1) + 2, 0, 0))] = cnt++;
143  }
144  }
145 
146  // Face 1
147  for (int i = 2; i <= P; ++i)
148  {
149  for (int j = 1; j <= R-i; ++j)
150  {
151  m_map[Mode(i, 0, j, GetTetMode(i, 0, j))] = cnt++;
152  }
153  }
154 
155  // Face 2
156  for (int i = 2; i <= Q; ++i)
157  {
158  for (int j = 1; j <= R-i; ++j)
159  {
160  m_map[Mode(1, i, j, GetTetMode(0, i, j))] = cnt++;
161  }
162  }
163 
164  // Face 3
165  for (int i = 2; i <= P; ++i)
166  {
167  for (int j = 1; j <= R-i; ++j)
168  {
169  m_map[Mode(i, 1, j, GetTetMode(i, 0, j))] = cnt++;
170  }
171  }
172 
173  // Face 4
174  for (int i = 2; i <= Q; ++i)
175  {
176  for (int j = 1; j <= R-i; ++j)
177  {
178  m_map[Mode(0, i, j, GetTetMode(0, i, j))] = cnt++;
179  }
180  }
181 
182  // Interior (tetrahedral modes)
183  for (int i = 2; i <= P+1; ++i)
184  {
185  for (int j = 1; j <= Q-i+1; ++j)
186  {
187  for (int k = 1; k <= R-i-j+1; ++k)
188  {
189  // need to go to j+1-th mode in the 'b' direction to
190  // select correct modified_a mode
191  m_map[Mode(i, j+1, k, GetTetMode(i-1, j, k))] = cnt++;
192  }
193  }
194  }
195 
196  ASSERTL0(m_map.size() == m_ncoeffs,
197  "Duplicate coefficient entries in map");
198 
200  for (it = m_map.begin(); it != m_map.end(); ++it)
201  {
202  const int p = it->first.get<0>();
203  const int q = it->first.get<1>();
204  const int r = it->first.get<2>();
205  const int rp = it->first.get<3>();
206  if (m_idxMap.find(p) == m_idxMap.end())
207  {
208  m_idxMap[p] = map<int, map<int, pair<int, int> > >();
209  }
210 
211  if (m_idxMap[p].find(q) == m_idxMap[p].end())
212  {
213  m_idxMap[p][q] = map<int, pair<int, int> >();
214  }
215 
216  if (m_idxMap[p][q].find(r) == m_idxMap[p][q].end())
217  {
218  m_idxMap[p][q][r] = pair<int, int>(it->second, rp);
219  }
220  }
221  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
boost::tuple< unsigned int, unsigned int, unsigned int, unsigned int > Mode
Definition: StdPyrExp.h:49
std::map< Mode, unsigned int, cmpop > m_map
Definition: StdPyrExp.h:255
int getNumberOfCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:232
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
int GetTetMode(int I, int J, int K)
Number tetrahedral modes in r-direction. Much the same as StdTetExp::GetTetMode but slightly simplifi...
Definition: StdPyrExp.cpp:1923
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:315
StdExpansion()
Default Constructor.
std::map< int, std::map< int, std::map< int, std::pair< int, int > > > > m_idxMap
Definition: StdPyrExp.h:256
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 223 of file StdPyrExp.cpp.

224  : StdExpansion (T),
225  StdExpansion3D(T)
226  {
227  }
StdExpansion()
Default Constructor.
Nektar::StdRegions::StdPyrExp::~StdPyrExp ( )

Definition at line 231 of file StdPyrExp.cpp.

232  {
233  }

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 1923 of file StdPyrExp.cpp.

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

Referenced by StdPyrExp().

1924  {
1925  const int R = m_base[2]->GetNumModes();
1926  int i, j, cnt = 0;
1927  for (i = 0; i < I; ++i)
1928  {
1929  cnt += (R-i)*(R-i+1)/2;
1930  }
1931 
1932  i = R-I;
1933  for (j = 0; j < J; ++j)
1934  {
1935  cnt += i;
1936  i--;
1937  }
1938 
1939  return cnt + K;
1940  }
Array< OneD, LibUtilities::BasisSharedPtr > m_base
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 369 of file StdPyrExp.cpp.

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

Referenced by v_FillMode().

371  {
372  if (m_base[0]->Collocation() &&
373  m_base[1]->Collocation() &&
374  m_base[2]->Collocation())
375  {
377  m_base[1]->GetNumPoints()*
378  m_base[2]->GetNumPoints(),
379  inarray, 1, outarray, 1);
380  }
381  else
382  {
383  StdPyrExp::v_BwdTrans_SumFac(inarray,outarray);
384  }
385  }
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:229
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdPyrExp.cpp:390
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
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 390 of file StdPyrExp.cpp.

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

Referenced by v_BwdTrans().

393  {
394  Array<OneD, NekDouble> wsp;
395  v_BwdTrans_SumFacKernel(m_base[0]->GetBdata(),
396  m_base[1]->GetBdata(),
397  m_base[2]->GetBdata(),
398  inarray,outarray,wsp,
399  true,true,true);
400  }
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)
Definition: StdPyrExp.cpp:402
Array< OneD, LibUtilities::BasisSharedPtr > m_base
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 402 of file StdPyrExp.cpp.

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

Referenced by v_BwdTrans_SumFac().

412  {
413  const int Qx = m_base[0]->GetNumPoints();
414  const int Qy = m_base[1]->GetNumPoints();
415  const int Qz = m_base[2]->GetNumPoints();
416 
417  const NekDouble *bx = base0.get();
418  const NekDouble *by = base1.get();
419  const NekDouble *bz = base2.get();
420 
421  // Need to count coeffs for storage...
422  map<int, map<int, map<int, pair<int, int> > > >::iterator it_p;
423  map<int, map<int, pair<int, int> > > ::iterator it_q;
424  map<int, pair<int, int> > ::iterator it_r;
425 
426  int pqcnt = 0;
427  for (it_p = m_idxMap.begin(); it_p != m_idxMap.end(); ++it_p)
428  {
429  for (it_q = it_p->second.begin(); it_q != it_p->second.end(); ++it_q)
430  {
431  pqcnt++;
432  }
433  }
434 
435  Array<OneD, NekDouble> fpq(pqcnt);
436  Array<OneD, NekDouble> fp (m_base[0]->GetNumModes());
437  int i ,j, k, s = 0, cnt = 0, cnt2 = 0;
438 
439  for (k = 0; k < Qz; ++k)
440  {
441  NekDouble bz1 = bz[k+Qz];
442 
443  cnt = 0;
444  for (it_p = m_idxMap.begin(); it_p != m_idxMap.end(); ++it_p)
445  {
446  for (it_q = it_p->second.begin(); it_q != it_p->second.end(); ++it_q)
447  {
448  NekDouble sum = 0.0;
449  for (it_r = it_q->second.begin(); it_r != it_q->second.end(); ++it_r)
450  {
451  sum += inarray[it_r->second.first] * bz[k + Qz*it_r->second.second];
452  }
453  fpq[cnt++] = sum;
454  }
455  }
456 
457  for (j = 0; j < Qy; ++j)
458  {
459  NekDouble by0 = bz1*by[j];
460  NekDouble by1 = bz1*by[j+Qy];
461 
462  cnt = cnt2 = 0;
463  for (it_p = m_idxMap.begin(); it_p != m_idxMap.end(); ++it_p)
464  {
465  NekDouble sum = 0.0;
466  for (it_q = it_p->second.begin(); it_q != it_p->second.end(); ++it_q)
467  {
468  sum += by[j + Qy*it_q->first] * fpq[cnt++];
469  }
470  fp[cnt2++] = sum;
471  }
472 
473  for (i = 0; i < Qx; ++i, ++s)
474  {
475  cnt2 = 0;
476  NekDouble sum = 0.0;
477  for (it_p = m_idxMap.begin(); it_p != m_idxMap.end(); ++it_p)
478  {
479  sum += bx[i + Qx*it_p->first] * fp[cnt2++];
480  }
481  sum += inarray[4]*(by1*(bx[i] + bx[i+Qx]) + by0*bx[i+Qx]);
482  outarray[s] = sum;
483  }
484  }
485  }
486  }
double NekDouble
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
Array< OneD, LibUtilities::BasisSharedPtr > m_base
std::map< int, std::map< int, std::map< int, std::pair< int, int > > > > m_idxMap
Definition: StdPyrExp.h:256
int Nektar::StdRegions::StdPyrExp::v_CalcNumberOfCoefficients ( const std::vector< unsigned int > &  nummodes,
int &  modes_offset 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1036 of file StdPyrExp.cpp.

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

1039  {
1041  nummodes[modes_offset],
1042  nummodes[modes_offset+1],
1043  nummodes[modes_offset+2]);
1044 
1045  modes_offset += 3;
1046  return nmodes;
1047  }
int getNumberOfCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:232
DNekMatSharedPtr Nektar::StdRegions::StdPyrExp::v_CreateStdMatrix ( const StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PyrExp.

Definition at line 1913 of file StdPyrExp.cpp.

References v_GenMatrix().

1914  {
1915  return v_GenMatrix(mkey);
1916  }
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey)
Definition: StdPyrExp.cpp:1908
const LibUtilities::BasisKey Nektar::StdRegions::StdPyrExp::v_DetFaceBasisKey ( const int  i,
const int  k 
) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 998 of file StdPyrExp.cpp.

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

1000  {
1001  ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
1002  ASSERTL2(k >= 0 && k <= 1, "basis key id is out of range");
1003 
1004  switch(i)
1005  {
1006  case 0:
1007  {
1008  return EvaluateQuadFaceBasisKey(k,
1009  m_base[k]->GetBasisType(),
1010  m_base[k]->GetNumPoints(),
1011  m_base[k]->GetNumModes());
1012 
1013  }
1014  case 1:
1015  case 3:
1016  {
1017  return EvaluateTriFaceBasisKey(k,
1018  m_base[2*k]->GetBasisType(),
1019  m_base[2*k]->GetNumPoints(),
1020  m_base[2*k]->GetNumModes());
1021  }
1022  case 2:
1023  case 4:
1024  {
1025  return EvaluateTriFaceBasisKey(k,
1026  m_base[k+1]->GetBasisType(),
1027  m_base[k+1]->GetNumPoints(),
1028  m_base[k+1]->GetNumModes());
1029  }
1030  }
1031 
1032  // Should never get here.
1034  }
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:229
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)
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:240
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.
Array< OneD, LibUtilities::BasisSharedPtr > m_base
LibUtilities::ShapeType Nektar::StdRegions::StdPyrExp::v_DetShapeType ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 891 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 864 of file StdPyrExp.cpp.

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

865  {
866  Array<OneD, NekDouble> tmp(m_ncoeffs, 0.0);
867  tmp[mode] = 1.0;
868  v_BwdTrans(tmp, outarray);
869  }
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Backward transformation is evaluated at the quadrature points.
Definition: StdPyrExp.cpp:369
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 501 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().

503  {
504  v_IProductWRTBase(inarray,outarray);
505 
506  // get Mass matrix inverse
507  StdMatrixKey masskey(eInvMass,DetShapeType(),*this);
508  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
509 
510  // copy inarray in case inarray == outarray
511  DNekVec in (m_ncoeffs, outarray);
512  DNekVec out(m_ncoeffs, outarray, eWrapper);
513 
514  out = (*matsys)*in;
515  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:470
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:700
NekVector< NekDouble > DNekVec
Definition: NekTypeDefs.hpp:49
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.
Definition: StdPyrExp.cpp:536
DNekMatSharedPtr Nektar::StdRegions::StdPyrExp::v_GenMatrix ( const StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PyrExp.

Definition at line 1908 of file StdPyrExp.cpp.

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

Referenced by v_CreateStdMatrix().

1909  {
1910  return CreateGeneralMatrix(mkey);
1911  }
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
void Nektar::StdRegions::StdPyrExp::v_GetBoundaryMap ( Array< OneD, unsigned int > &  outarray)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1884 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().

1885  {
1888  "BasisType is not a boundary interior form");
1891  "BasisType is not a boundary interior form");
1894  "BasisType is not a boundary interior form");
1895 
1896  int idx = 0, nBndry = v_NumBndryCoeffs();
1897 
1898  for (idx = 0; idx < nBndry; ++idx)
1899  {
1900  maparray[idx] = idx;
1901  }
1902  }
Principle Modified Functions .
Definition: BasisType.h:51
Principle Modified Functions .
Definition: BasisType.h:49
virtual int v_NumBndryCoeffs() const
Definition: StdPyrExp.cpp:896
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
Lagrange for SEM basis .
Definition: BasisType.h:53
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
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 836 of file StdPyrExp.cpp.

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

839  {
840  Array<OneD, const NekDouble> etaBar_x = m_base[0]->GetZ();
841  Array<OneD, const NekDouble> eta_y = m_base[1]->GetZ();
842  Array<OneD, const NekDouble> eta_z = m_base[2]->GetZ();
843  int Qx = GetNumPoints(0);
844  int Qy = GetNumPoints(1);
845  int Qz = GetNumPoints(2);
846 
847  // Convert collapsed coordinates into cartesian coordinates: eta --> xi
848  for (int k = 0; k < Qz; ++k )
849  {
850  for (int j = 0; j < Qy; ++j)
851  {
852  for (int i = 0; i < Qx; ++i)
853  {
854  int s = i + Qx*(j + Qy*k);
855 
856  xi_z[s] = eta_z[k];
857  xi_y[s] = (1.0 + eta_y[j]) * (1.0 - eta_z[k]) / 2.0 - 1.0;
858  xi_x[s] = (1.0 + etaBar_x[i]) * (1.0 - eta_z[k]) / 2.0 - 1.0;
859  }
860  }
861  }
862  }
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:229
Array< OneD, LibUtilities::BasisSharedPtr > m_base
LibUtilities::BasisType Nektar::StdRegions::StdPyrExp::v_GetEdgeBasisType ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1049 of file StdPyrExp.cpp.

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

1050  {
1051  ASSERTL2(i >= 0 && i <= 7, "edge id is out of range");
1052  if (i == 0 || i == 2)
1053  {
1054  return GetBasisType(0);
1055  }
1056  else if (i == 1 || i == 3)
1057  {
1058  return GetBasisType(1);
1059  }
1060  else
1061  {
1062  return GetBasisType(2);
1063  }
1064  }
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:240
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 1614 of file StdPyrExp.cpp.

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

1619  {
1620  int i;
1621  bool signChange;
1622  const int P = m_base[0]->GetNumModes() - 2;
1623  const int Q = m_base[1]->GetNumModes() - 2;
1624  const int R = m_base[2]->GetNumModes() - 2;
1625  const int nEdgeIntCoeffs = v_GetEdgeNcoeffs(eid) - 2;
1626 
1627  if (maparray.num_elements() != nEdgeIntCoeffs)
1628  {
1629  maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
1630  }
1631 
1632  if(signarray.num_elements() != nEdgeIntCoeffs)
1633  {
1634  signarray = Array<OneD, int>(nEdgeIntCoeffs,1);
1635  }
1636  else
1637  {
1638  fill(signarray.get(), signarray.get()+nEdgeIntCoeffs, 1);
1639  }
1640 
1641  // If edge is oriented backwards, change sign of modes which have
1642  // degree 2n+1, n >= 1.
1643  signChange = edgeOrient == eBackwards;
1644 
1645  int offset = 5;
1646 
1647  switch (eid)
1648  {
1649  case 0:
1650  break;
1651  case 1:
1652  offset += P;
1653  break;
1654  case 2:
1655  offset += P+Q;
1656  break;
1657  case 3:
1658  offset += 2*P+Q;
1659  break;
1660  case 4:
1661  offset += 2*(P+Q);
1662  break;
1663  case 5:
1664  offset += 2*(P+Q)+R;
1665  break;
1666  case 6:
1667  offset += 2*(P+Q+R);
1668  break;
1669  case 7:
1670  offset += 2*(P+Q)+3*R;
1671  break;
1672  default:
1673  ASSERTL0(false, "Edge not defined.");
1674  break;
1675  }
1676 
1677  for (i = 0; i < nEdgeIntCoeffs; ++i)
1678  {
1679  maparray[i] = offset + i;
1680  }
1681 
1682  if (signChange)
1683  {
1684  for (i = 1; i < nEdgeIntCoeffs; i += 2)
1685  {
1686  signarray[i] = -1;
1687  }
1688  }
1689  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
virtual int v_GetEdgeNcoeffs(const int i) const
Definition: StdPyrExp.cpp:916
Array< OneD, LibUtilities::BasisSharedPtr > m_base
int Nektar::StdRegions::StdPyrExp::v_GetEdgeNcoeffs ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 916 of file StdPyrExp.cpp.

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

Referenced by v_GetEdgeInteriorMap().

917  {
918  ASSERTL2(i >= 0 && i <= 7, "edge id is out of range");
919 
920  if (i == 0 || i == 2)
921  {
922  return GetBasisNumModes(0);
923  }
924  else if (i == 1 || i == 3)
925  {
926  return GetBasisNumModes(1);
927  }
928  else
929  {
930  return GetBasisNumModes(2);
931  }
932  }
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:178
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:240
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 1691 of file StdPyrExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 954 of file StdPyrExp.cpp.

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

Referenced by v_GetFaceInteriorMap(), and v_GetFaceToElementMap().

955  {
956  ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
957 
958  int P = m_base[0]->GetNumModes()-1;
959  int Q = m_base[1]->GetNumModes()-1;
960  int R = m_base[2]->GetNumModes()-1;
961 
962  if (i == 0)
963  {
964  return (P-1)*(Q-1);
965  }
966  else if (i == 1 || i == 3)
967  {
968  return (P-1) * (2*(R-1) - (P-1) - 1) / 2;
969  }
970  else
971  {
972  return (Q-1) * (2*(R-1) - (Q-1) - 1) / 2;
973  }
974  }
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:240
Array< OneD, LibUtilities::BasisSharedPtr > m_base
int Nektar::StdRegions::StdPyrExp::v_GetFaceNcoeffs ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 934 of file StdPyrExp.cpp.

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

935  {
936  ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
937 
938  if (i == 0)
939  {
940  return GetBasisNumModes(0)*GetBasisNumModes(1);
941  }
942  else if (i == 1 || i == 3)
943  {
944  int P = GetBasisNumModes(0)-1, Q = GetBasisNumModes(2)-1;
945  return Q+1 + (P*(1 + 2*Q - P))/2;
946  }
947  else
948  {
949  int P = GetBasisNumModes(1)-1, Q = GetBasisNumModes(2)-1;
950  return Q+1 + (P*(1 + 2*Q - P))/2;
951  }
952  }
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:178
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:240
int Nektar::StdRegions::StdPyrExp::v_GetFaceNumPoints ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 976 of file StdPyrExp.cpp.

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

977  {
978  ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
979 
980  if (i == 0)
981  {
982  return m_base[0]->GetNumPoints()*
983  m_base[1]->GetNumPoints();
984  }
985  else if (i == 1 || i == 3)
986  {
987  return m_base[0]->GetNumPoints()*
988  m_base[2]->GetNumPoints();
989  }
990  else
991  {
992  return m_base[1]->GetNumPoints()*
993  m_base[2]->GetNumPoints();
994  }
995  }
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:240
Array< OneD, LibUtilities::BasisSharedPtr > m_base
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 1071 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().

1078  {
1080  "Method only implemented if BasisType is identical"
1081  "in x and y directions");
1084  "Method only implemented for Modified_A BasisType"
1085  "(x and y direction) and Modified_C BasisType (z "
1086  "direction)");
1087 
1088  int i, j, k, p, q, r, nFaceCoeffs;
1089  int nummodesA, nummodesB;
1090 
1091  int order0 = m_base[0]->GetNumModes();
1092  int order1 = m_base[1]->GetNumModes();
1093  int order2 = m_base[2]->GetNumModes();
1094 
1095  switch (fid)
1096  {
1097  case 0:
1098  nummodesA = order0;
1099  nummodesB = order1;
1100  break;
1101  case 1:
1102  case 3:
1103  nummodesA = order0;
1104  nummodesB = order2;
1105  break;
1106  case 2:
1107  case 4:
1108  nummodesA = order1;
1109  nummodesB = order2;
1110  break;
1111  }
1112 
1113  bool CheckForZeroedModes = false;
1114 
1115  if (P == -1)
1116  {
1117  P = nummodesA;
1118  Q = nummodesB;
1119  nFaceCoeffs = GetFaceNcoeffs(fid);
1120  }
1121  else if (fid > 0)
1122  {
1123  nFaceCoeffs = P*(2*Q-P+1)/2;
1124  CheckForZeroedModes = true;
1125  }
1126  else
1127  {
1128  nFaceCoeffs = P*Q;
1129  CheckForZeroedModes = true;
1130  }
1131 
1132  // Allocate the map array and sign array; set sign array to ones (+)
1133  if (maparray.num_elements() != nFaceCoeffs)
1134  {
1135  maparray = Array<OneD, unsigned int>(nFaceCoeffs);
1136  }
1137 
1138  if (signarray.num_elements() != nFaceCoeffs)
1139  {
1140  signarray = Array<OneD, int>(nFaceCoeffs,1);
1141  }
1142  else
1143  {
1144  fill(signarray.get(), signarray.get() + nFaceCoeffs, 1);
1145  }
1146 
1147  // Set up an array indexing for quads, since the ordering may need
1148  // to be transposed.
1149  Array<OneD, int> arrayindx(nFaceCoeffs,-1);
1150 
1151  if (fid == 0)
1152  {
1153  for (i = 0; i < Q; i++)
1154  {
1155  for (j = 0; j < P; j++)
1156  {
1157  if (faceOrient < 9)
1158  {
1159  arrayindx[i*P+j] = i*P+j;
1160  }
1161  else
1162  {
1163  arrayindx[i*P+j] = j*Q+i;
1164  }
1165  }
1166  }
1167  }
1168 
1169  // Set up ordering inside each 2D face. Also for triangular faces,
1170  // populate signarray.
1171  int cnt = 0, cnt2;
1172  switch (fid)
1173  {
1174  case 0: // Bottom quad
1175 
1176  // Fill in vertices
1177  maparray[arrayindx[0]] = 0;
1178  maparray[arrayindx[1]] = 1;
1179  maparray[arrayindx[P+1]] = 2;
1180  maparray[arrayindx[P]] = 3;
1181 
1182  // Edge 0
1183  cnt = 5;
1184  for (p = 2; p < P; ++p)
1185  {
1186  maparray[arrayindx[p]] = p-2 + cnt;
1187  }
1188 
1189  // Edge 1
1190  cnt += P-2;
1191  for (q = 2; q < Q; ++q)
1192  {
1193  maparray[arrayindx[q*P+1]] = q-2 + cnt;
1194  }
1195 
1196  // Edge 2
1197  cnt += Q-2;
1198  for (p = 2; p < P; ++p)
1199  {
1200  maparray[arrayindx[P+p]] = p-2 + cnt;
1201  }
1202 
1203  // Edge 3
1204  cnt += P-2;
1205  for (q = 2; q < Q; ++q)
1206  {
1207  maparray[arrayindx[q*P]] = q-2 + cnt;
1208  }
1209 
1210  // Interior
1211  cnt += Q-2 + 4*(P-2);
1212  for (q = 2; q < Q; ++q)
1213  {
1214  for (p = 2; p < P; ++p)
1215  {
1216  maparray[arrayindx[q*P+p]] = cnt + (q-2)*P+(p-2);
1217  }
1218  }
1219 
1220 
1221  break;
1222 
1223  case 1: // Left triangle
1224  // Vertices
1225  maparray[0] = 0;
1226  maparray[1] = 4;
1227  maparray[Q] = 1;
1228 
1229  // Edge 0 (pyramid edge 0)
1230  cnt = 5;
1231  q = 2*Q-1;
1232  for (p = 2; p < P; q += Q-p, ++p)
1233  {
1234  maparray[q] = cnt++;
1235  if ((int)faceOrient == 7)
1236  {
1237  signarray[q] = p % 2 ? -1 : 1;
1238  }
1239  }
1240 
1241  // Edge 1 (pyramid edge 5)
1242  cnt = 5 + 2*(order0-2) + 2*(order1-2) + (order2-2);
1243  for (q = 2; q < Q; ++q)
1244  {
1245  maparray[q] = cnt++;
1246  if ((int)faceOrient == 7)
1247  {
1248  signarray[q] = q % 2 ? -1 : 1;
1249  }
1250  }
1251 
1252  // Edge 2 (pyramid edge 4)
1253  cnt = 5 + 2*(order0-2) + 2*(order1-2);
1254  for (q = 2; q < Q; ++q)
1255  {
1256  maparray[Q+q-1] = cnt++;
1257  if ((int)faceOrient == 7)
1258  {
1259  signarray[Q+q-1] = q % 2 ? -1 : 1;
1260  }
1261  }
1262 
1263  // Interior
1264  cnt = 5 + 2*(order0-2) + 2*(order1-2) + 4*(order2-2)
1265  + v_GetFaceIntNcoeffs(0);
1266  cnt2 = 2*Q + 1;
1267  for (p = 2; p < P; ++p)
1268  {
1269  for (r = 2; r < Q-p; ++r)
1270  {
1271  maparray[cnt2] = cnt++;
1272  if ((int)faceOrient == 7 && p > 1)
1273  {
1274  signarray[cnt2++] = p % 2 ? -1 : 1;
1275  }
1276  }
1277  cnt2++;
1278  }
1279  break;
1280 
1281  case 2:
1282  // Vertices
1283  maparray[0] = 1;
1284  maparray[1] = 4;
1285  maparray[Q] = 2;
1286 
1287  // Edge 0 (pyramid edge 1)
1288  cnt = 5 + (order0-2);
1289  q = 2*Q-1;
1290  for (p = 2; p < P; q += Q-p, ++p)
1291  {
1292  maparray[q] = cnt++;
1293  if ((int)faceOrient == 7)
1294  {
1295  signarray[q] = p % 2 ? -1 : 1;
1296  }
1297  }
1298 
1299  // Edge 1 (pyramid edge 6)
1300  cnt = 5 + 2*(order0-2) + 2*(order1-2) + 2*(order2-2);
1301  for (q = 2; q < Q; ++q)
1302  {
1303  maparray[q] = cnt++;
1304  if ((int)faceOrient == 7)
1305  {
1306  signarray[q] = q % 2 ? -1 : 1;
1307  }
1308  }
1309 
1310  // Edge 2 (pyramid edge 5)
1311  cnt = 5 + 2*(order0-2) + 2*(order1-2) + (order2-2);
1312  for (q = 2; q < Q; ++q)
1313  {
1314  maparray[Q+q-1] = cnt++;
1315  if ((int)faceOrient == 7)
1316  {
1317  signarray[Q+q-1] = q % 2 ? -1 : 1;
1318  }
1319  }
1320 
1321  // Interior
1322  cnt = 5 + 2*(order0-2) + 2*(order1-2) + 4*(order2-2)
1324  cnt2 = 2*Q + 1;
1325  for (p = 2; p < P; ++p)
1326  {
1327  for (r = 2; r < Q-p; ++r)
1328  {
1329  maparray[cnt2] = cnt++;
1330  if ((int)faceOrient == 7 && p > 1)
1331  {
1332  signarray[cnt2++] = p % 2 ? -1 : 1;
1333  }
1334  }
1335  cnt2++;
1336  }
1337  break;
1338 
1339  case 3: // Right triangle
1340  // Vertices
1341  maparray[0] = 3;
1342  maparray[1] = 4;
1343  maparray[Q] = 2;
1344 
1345  // Edge 0 (pyramid edge 2)
1346  cnt = 5 + (order0-2) + (order1-2);
1347  q = 2*Q-1;
1348  for (p = 2; p < P; q += Q-p, ++p)
1349  {
1350  maparray[q] = cnt++;
1351  if ((int)faceOrient == 7)
1352  {
1353  signarray[q] = p % 2 ? -1 : 1;
1354  }
1355  }
1356 
1357  // Edge 1 (pyramid edge 6)
1358  cnt = 5 + 2*(order0-2) + 2*(order1-2) + 2*(order2-2);
1359  for (q = 2; q < Q; ++q)
1360  {
1361  maparray[q] = cnt++;
1362  if ((int)faceOrient == 7)
1363  {
1364  signarray[q] = q % 2 ? -1 : 1;
1365  }
1366  }
1367 
1368  // Edge 2 (pyramid edge 7)
1369  cnt = 5 + 2*(order0-2) + 2*(order1-2) + 3*(order2-2);
1370  for (q = 2; q < Q; ++q)
1371  {
1372  maparray[Q+q-1] = cnt++;
1373  if ((int)faceOrient == 7)
1374  {
1375  signarray[Q+q-1] = q % 2 ? -1 : 1;
1376  }
1377  }
1378 
1379  // Interior
1380  cnt = 5 + 2*(order0-2) + 2*(order1-2) + 4*(order2-2)
1382  + v_GetFaceIntNcoeffs(2);
1383  cnt2 = 2*Q + 1;
1384  for (p = 2; p < P; ++p)
1385  {
1386  for (r = 2; r < Q-p; ++r)
1387  {
1388  maparray[cnt2] = cnt++;
1389  if ((int)faceOrient == 7 && p > 1)
1390  {
1391  signarray[cnt2++] = p % 2 ? -1 : 1;
1392  }
1393  }
1394  cnt2++;
1395  }
1396  break;
1397 
1398  case 4: // Rear tri
1399  // Vertices
1400  maparray[0] = 0;
1401  maparray[1] = 4;
1402  maparray[Q] = 3;
1403 
1404  // Edge 0 (pyramid edge 3)
1405  cnt = 5 + 2*(order0-2) + (order1-2);
1406  q = 2*Q-1;
1407  for (p = 2; p < P; q += Q-p, ++p)
1408  {
1409  maparray[q] = cnt++;
1410  if ((int)faceOrient == 7)
1411  {
1412  signarray[q] = p % 2 ? -1 : 1;
1413  }
1414  }
1415 
1416  // Edge 1 (pyramid edge 7)
1417  cnt = 5 + 2*(order0-2) + 2*(order1-2) + 3*(order2-2);
1418  for (q = 2; q < Q; ++q)
1419  {
1420  maparray[q] = cnt++;
1421  if ((int)faceOrient == 7)
1422  {
1423  signarray[q] = q % 2 ? -1 : 1;
1424  }
1425  }
1426 
1427  // Edge 2 (pyramid edge 4)
1428  cnt = 5 + 2*(order0-2) + 2*(order1-2);
1429  for (q = 2; q < Q; ++q)
1430  {
1431  maparray[Q+q-1] = cnt++;
1432  if ((int)faceOrient == 7)
1433  {
1434  signarray[Q+q-1] = q % 2 ? -1 : 1;
1435  }
1436  }
1437 
1438  // Interior
1439  cnt = 5 + 2*(order0-2) + 2*(order1-2) + 4*(order2-2)
1442  cnt2 = 2*Q + 1;
1443  for (p = 2; p < P; ++p)
1444  {
1445  for (r = 2; r < Q-p; ++r)
1446  {
1447  maparray[cnt2] = cnt++;
1448  if ((int)faceOrient == 7 && p > 1)
1449  {
1450  signarray[cnt2++] = p % 2 ? -1 : 1;
1451  }
1452  }
1453  cnt2++;
1454  }
1455  break;
1456 
1457  default:
1458  ASSERTL0(false, "Face to element map unavailable.");
1459  }
1460 
1461  if (fid > 0)
1462  {
1463 
1464  if(CheckForZeroedModes)
1465  {
1466  // zero signmap and set maparray to zero if elemental
1467  // modes are not as large as face modesl
1468  int idx = 0;
1469  for (j = 0; j < P; ++j)
1470  {
1471  idx += Q-j;
1472  for (k = Q-j; k < Q-j; ++k)
1473  {
1474  signarray[idx] = 0.0;
1475  maparray[idx++] = maparray[0];
1476  }
1477  }
1478 
1479  for (j = P; j < P; ++j)
1480  {
1481  for (k = 0; k < Q-j; ++k)
1482  {
1483  signarray[idx] = 0.0;
1484  maparray[idx++] = maparray[0];
1485  }
1486  }
1487  }
1488 
1489  // Triangles only have one possible orientation (base
1490  // direction reversed); swap edge modes.
1491  if ((int)faceOrient == 7)
1492  {
1493  swap(maparray[0], maparray[P]);
1494  for (i = 1; i < P-1; ++i)
1495  {
1496  swap(maparray[i+1], maparray[P+i]);
1497  }
1498  }
1499  }
1500  else
1501  {
1502  if(CheckForZeroedModes)
1503  {
1504  // zero signmap and set maparray to zero if elemental
1505  // modes are not as large as face modesl
1506  for (j = 0; j < P; ++j)
1507  {
1508  for (k = Q; k < Q; ++k)
1509  {
1510  signarray[arrayindx[j+k*P]] = 0.0;
1511  maparray[arrayindx[j+k*P]] = maparray[0];
1512  }
1513  }
1514 
1515  for (j = P; j < P; ++j)
1516  {
1517  for (k = 0; k < Q; ++k)
1518  {
1519  signarray[arrayindx[j+k*P]] = 0.0;
1520  maparray[arrayindx[j+k*P]] = maparray[0];
1521  }
1522  }
1523  }
1524 
1525  // The code below is exactly the same as that taken from
1526  // StdHexExp and reverses the 'b' and 'a' directions as
1527  // appropriate (1st and 2nd if statements respectively) in
1528  // quadrilateral faces.
1529  if (faceOrient == 6 || faceOrient == 8 ||
1530  faceOrient == 11 || faceOrient == 12)
1531  {
1532  if (faceOrient < 9)
1533  {
1534  for (i = 3; i < Q; i += 2)
1535  {
1536  for (j = 0; j < P; j++)
1537  {
1538  signarray[arrayindx[i*P+j]] *= -1;
1539  }
1540  }
1541 
1542  for (i = 0; i < P; i++)
1543  {
1544  swap(maparray [i], maparray [i+P]);
1545  swap(signarray[i], signarray[i+P]);
1546  }
1547  }
1548  else
1549  {
1550  for (i = 0; i < Q; i++)
1551  {
1552  for (j = 3; j < P; j += 2)
1553  {
1554  signarray[arrayindx[i*P+j]] *= -1;
1555  }
1556  }
1557 
1558  for (i = 0; i < Q; i++)
1559  {
1560  swap (maparray [i], maparray [i+Q]);
1561  swap (signarray[i], signarray[i+Q]);
1562  }
1563  }
1564  }
1565 
1566  if (faceOrient == 7 || faceOrient == 8 ||
1567  faceOrient == 10 || faceOrient == 12)
1568  {
1569  if (faceOrient < 9)
1570  {
1571  for (i = 0; i < Q; i++)
1572  {
1573  for (j = 3; j < P; j += 2)
1574  {
1575  signarray[arrayindx[i*P+j]] *= -1;
1576  }
1577  }
1578 
1579  for(i = 0; i < Q; i++)
1580  {
1581  swap(maparray [i*P], maparray [i*P+1]);
1582  swap(signarray[i*P], signarray[i*P+1]);
1583  }
1584  }
1585  else
1586  {
1587  for (i = 3; i < Q; i += 2)
1588  {
1589  for (j = 0; j < P; j++)
1590  {
1591  signarray[arrayindx[i*P+j]] *= -1;
1592  }
1593  }
1594 
1595  for (i = 0; i < P; i++)
1596  {
1597  swap(maparray [i*Q], maparray [i*Q+1]);
1598  swap(signarray[i*Q], signarray[i*Q+1]);
1599  }
1600  }
1601  }
1602  }
1603  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
Principle Modified Functions .
Definition: BasisType.h:51
Principle Modified Functions .
Definition: BasisType.h:49
LibUtilities::BasisType GetEdgeBasisType(const int i) const
This function returns the type of expansion basis on the i-th edge.
Definition: StdExpansion.h:413
virtual int v_GetFaceIntNcoeffs(const int i) const
Definition: StdPyrExp.cpp:954
int GetFaceNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th face.
Definition: StdExpansion.h:354
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdPyrExp::v_GetInteriorMap ( Array< OneD, unsigned int > &  outarray)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1856 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().

1857  {
1860  "BasisType is not a boundary interior form");
1863  "BasisType is not a boundary interior form");
1866  "BasisType is not a boundary interior form");
1867 
1868  const int nBndCoeffs = v_NumBndryCoeffs();
1869  const int nIntCoeffs = m_ncoeffs - NumBndryCoeffs();
1870 
1871  if (outarray.num_elements() != nIntCoeffs)
1872  {
1873  outarray = Array<OneD, unsigned int>(nIntCoeffs);
1874  }
1875 
1876  // Loop over all interior modes.
1877  int p, idx = 0;
1878  for (p = nBndCoeffs; p < m_ncoeffs; ++p)
1879  {
1880  outarray[idx++] = p;
1881  }
1882  }
Principle Modified Functions .
Definition: BasisType.h:51
Principle Modified Functions .
Definition: BasisType.h:49
virtual int v_NumBndryCoeffs() const
Definition: StdPyrExp.cpp:896
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
Lagrange for SEM basis .
Definition: BasisType.h:53
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
int Nektar::StdRegions::StdPyrExp::v_GetNedges ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 881 of file StdPyrExp.cpp.

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 886 of file StdPyrExp.cpp.

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 876 of file StdPyrExp.cpp.

877  {
878  return 5;
879  }
int Nektar::StdRegions::StdPyrExp::v_GetVertexMap ( int  localVertexId,
bool  useCoeffPacking = false 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1605 of file StdPyrExp.cpp.

References ASSERTL0, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_C, and Nektar::StdRegions::StdExpansion::GetEdgeBasisType().

1606  {
1610  "Mapping not defined for this type of basis");
1611  return vId;
1612  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
Principle Modified Functions .
Definition: BasisType.h:51
Principle Modified Functions .
Definition: BasisType.h:49
LibUtilities::BasisType GetEdgeBasisType(const int i) const
This function returns the type of expansion basis on the i-th edge.
Definition: StdExpansion.h:413
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 536 of file StdPyrExp.cpp.

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

Referenced by v_FwdTrans().

539  {
540  if (m_base[0]->Collocation() &&
541  m_base[1]->Collocation() &&
542  m_base[2]->Collocation())
543  {
544  MultiplyByStdQuadratureMetric(inarray, outarray);
545  }
546  else
547  {
548  StdPyrExp::v_IProductWRTBase_SumFac(inarray,outarray);
549  }
550  }
void MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:949
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
Definition: StdPyrExp.cpp:552
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdPyrExp::v_IProductWRTBase_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
bool  multiplybyweights = true 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 552 of file StdPyrExp.cpp.

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

Referenced by v_IProductWRTBase().

556  {
557  Array<OneD, NekDouble> wsp;
558 
559  if(multiplybyweights)
560  {
561  Array<OneD, NekDouble> tmp(inarray.num_elements());
562 
563  v_MultiplyByStdQuadratureMetric(inarray, tmp);
564 
565  v_IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
566  m_base[1]->GetBdata(),
567  m_base[2]->GetBdata(),
568  tmp,outarray,wsp,
569  true,true,true);
570  }
571  else
572  {
573  v_IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
574  m_base[1]->GetBdata(),
575  m_base[2]->GetBdata(),
576  inarray,outarray,wsp,
577  true,true,true);
578  }
579  }
void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdPyrExp.cpp:1942
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)
Definition: StdPyrExp.cpp:581
Array< OneD, LibUtilities::BasisSharedPtr > m_base
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 581 of file StdPyrExp.cpp.

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

Referenced by v_IProductWRTBase_SumFac().

591  {
592  int i, j, k, s;
593  int Qx = m_base[0]->GetNumPoints();
594  int Qy = m_base[1]->GetNumPoints();
595  int Qz = m_base[2]->GetNumPoints();
596 
597  const NekDouble *bx = base0.get();
598  const NekDouble *by = base1.get();
599  const NekDouble *bz = base2.get();
600 
601  map<int, map<int, map<int, pair<int, int> > > >::iterator it_p;
602  map<int, map<int, pair<int, int> > > ::iterator it_q;
603  map<int, pair<int, int> > ::iterator it_r;
604 
605  Array<OneD, NekDouble> f (Qy*Qz);
606  Array<OneD, NekDouble> fb(Qz);
607 
608  for (it_p = m_idxMap.begin(); it_p != m_idxMap.end(); ++it_p)
609  {
610  const int p = it_p->first;
611  s = 0;
612  for (k = 0; k < Qz; ++k)
613  {
614  for (j = 0; j < Qy; ++j)
615  {
616  NekDouble sum = 0.0;
617  for (i = 0; i < Qx; ++i, ++s)
618  {
619  sum += bx[i + Qx*p]*inarray[s];
620  }
621  f[j+Qy*k] = sum;
622  }
623  }
624 
625  for (it_q = it_p->second.begin(); it_q != it_p->second.end(); ++it_q)
626  {
627  const int q = it_q->first;
628 
629  for (k = 0; k < Qz; ++k)
630  {
631  NekDouble sum = 0.0;
632  for (j = 0; j < Qy; ++j)
633  {
634  sum += by[j + Qy*q]*f[j+Qy*k];
635  }
636  fb[k] = sum;
637  }
638 
639  for (it_r = it_q->second.begin(); it_r != it_q->second.end(); ++it_r)
640  {
641  const int rpqr = it_r->second.second;
642  NekDouble sum = 0.0;
643  for (k = 0; k < Qz; ++k)
644  {
645  sum += bz[k + Qz*rpqr]*fb[k];
646  }
647 
648  outarray[it_r->second.first] = sum;
649  }
650  }
651  }
652 
653  // Correct for top mode
654  s = 0;
655  for (k = 0; k < Qz; ++k)
656  {
657  for (j = 0; j < Qy; ++j)
658  {
659  for (i = 0; i < Qx; ++i, ++s)
660  {
661  outarray[4] += inarray[s] * bz[k+Qz]*(
662  bx[i+Qx]*by[j+Qy] +
663  bx[i+Qx]*by[j ] +
664  bx[i ]*by[j+Qy]);
665  }
666  }
667  }
668  }
double NekDouble
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
Array< OneD, LibUtilities::BasisSharedPtr > m_base
std::map< int, std::map< int, std::map< int, std::pair< int, int > > > > m_idxMap
Definition: StdPyrExp.h:256
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 670 of file StdPyrExp.cpp.

References v_IProductWRTDerivBase_SumFac().

674  {
675  StdPyrExp::v_IProductWRTDerivBase_SumFac(dir,inarray,outarray);
676  }
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdPyrExp.cpp:684
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 684 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().

688  {
689  int i;
690  int nquad0 = m_base[0]->GetNumPoints();
691  int nquad1 = m_base[1]->GetNumPoints();
692  int nquad2 = m_base[2]->GetNumPoints();
693  int nqtot = nquad0*nquad1*nquad2;
694 
695  Array<OneD, NekDouble> gfac0(nquad0);
696  Array<OneD, NekDouble> gfac1(nquad1);
697  Array<OneD, NekDouble> gfac2(nquad2);
698  Array<OneD, NekDouble> tmp0 (nqtot);
699  Array<OneD, NekDouble> wsp;
700 
701  const Array<OneD, const NekDouble>& z0 = m_base[0]->GetZ();
702  const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
703  const Array<OneD, const NekDouble>& z2 = m_base[2]->GetZ();
704 
705  // set up geometric factor: (1+z0)/2
706  for(i = 0; i < nquad0; ++i)
707  {
708  gfac0[i] = 0.5*(1+z0[i]);
709  }
710 
711  // set up geometric factor: (1+z1)/2
712  for(i = 0; i < nquad1; ++i)
713  {
714  gfac1[i] = 0.5*(1+z1[i]);
715  }
716 
717  // Set up geometric factor: 2/(1-z2)
718  for(i = 0; i < nquad2; ++i)
719  {
720  gfac2[i] = 2.0/(1-z2[i]);
721  }
722 
723  // Derivative in first/second direction is always scaled as follows
724  const int nq01 = nquad0*nquad1;
725  for(i = 0; i < nquad2; ++i)
726  {
727  Vmath::Smul(nq01, gfac2[i],
728  &inarray[0] + i*nq01, 1,
729  &tmp0 [0] + i*nq01, 1);
730  }
731 
732  MultiplyByStdQuadratureMetric(tmp0, tmp0);
733 
734  switch(dir)
735  {
736  case 0:
737  {
738  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
739  m_base[1]->GetBdata (),
740  m_base[2]->GetBdata (),
741  tmp0, outarray, wsp,
742  false, true, true);
743  break;
744  }
745  case 1:
746  {
747  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata (),
748  m_base[1]->GetDbdata(),
749  m_base[2]->GetBdata (),
750  tmp0, outarray, wsp,
751  true, false, true);
752  break;
753  }
754  case 2:
755  {
756  Array<OneD, NekDouble> tmp3(m_ncoeffs);
757  Array<OneD, NekDouble> tmp4(m_ncoeffs);
758 
759  // Scale eta_1 derivative by gfac0
760  for(i = 0; i < nquad1*nquad2; ++i)
761  {
762  Vmath::Vmul(nquad0, tmp0 .get() + i*nquad0, 1,
763  gfac0.get(), 1,
764  tmp0 .get() + i*nquad0, 1);
765  }
766  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
767  m_base[1]->GetBdata(),
768  m_base[2]->GetBdata(),
769  tmp0, tmp3, wsp,
770  false, true, true);
771 
772  // Scale eta_2 derivative by gfac1*gfac2
773  for(i = 0; i < nquad2; ++i)
774  {
775  Vmath::Smul(nq01, gfac2[i],
776  &inarray[0] + i*nq01, 1,
777  &tmp0 [0] + i*nq01, 1);
778  }
779  for(i = 0; i < nquad1*nquad2; ++i)
780  {
781  Vmath::Smul(nquad0, gfac1[i%nquad1],
782  &tmp0[0] + i*nquad0, 1,
783  &tmp0[0] + i*nquad0, 1);
784  }
785 
786  MultiplyByStdQuadratureMetric(tmp0, tmp0);
787  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
788  m_base[1]->GetDbdata(),
789  m_base[2]->GetBdata(),
790  tmp0, tmp4, wsp,
791  true, false, true);
792 
793  MultiplyByStdQuadratureMetric(inarray,tmp0);
794  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
795  m_base[1]->GetBdata(),
796  m_base[2]->GetDbdata(),
797  tmp0,outarray,wsp,
798  true, true, false);
799 
800  Vmath::Vadd(m_ncoeffs,&tmp3[0],1,&outarray[0],1,&outarray[0],1);
801  Vmath::Vadd(m_ncoeffs,&tmp4[0],1,&outarray[0],1,&outarray[0],1);
802  break;
803  }
804  default:
805  {
806  ASSERTL1(false, "input dir is out of range");
807  break;
808  }
809  }
810  }
void MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:949
void IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:285
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:169
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 816 of file StdPyrExp.cpp.

References Nektar::NekConstants::kNekZeroTol.

819  {
820  if (fabs(xi[2]-1.0) < NekConstants::kNekZeroTol)
821  {
822  // Very top point of the pyramid
823  eta[0] = -1.0;
824  eta[1] = -1.0;
825  eta[2] = xi[2];
826  }
827  else
828  {
829  // Below the line-singularity -- Common case
830  eta[2] = xi[2]; // eta_z = xi_z
831  eta[1] = 2.0*(1.0 + xi[1])/(1.0 - xi[2]) - 1.0;
832  eta[0] = 2.0*(1.0 + xi[0])/(1.0 - xi[2]) - 1.0;
833  }
834  }
static const NekDouble kNekZeroTol
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 1942 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().

1945  {
1946  int i, j;
1947 
1948  int nquad0 = m_base[0]->GetNumPoints();
1949  int nquad1 = m_base[1]->GetNumPoints();
1950  int nquad2 = m_base[2]->GetNumPoints();
1951 
1952  const Array<OneD, const NekDouble>& w0 = m_base[0]->GetW();
1953  const Array<OneD, const NekDouble>& w1 = m_base[1]->GetW();
1954  const Array<OneD, const NekDouble>& w2 = m_base[2]->GetW();
1955 
1956  const Array<OneD, const NekDouble>& z2 = m_base[2]->GetZ();
1957 
1958  // Multiply by integration constants in x-direction
1959  for(i = 0; i < nquad1*nquad2; ++i)
1960  {
1961  Vmath::Vmul(nquad0, inarray.get()+i*nquad0, 1,
1962  w0.get(), 1, outarray.get()+i*nquad0,1);
1963  }
1964 
1965  // Multiply by integration constants in y-direction
1966  for(j = 0; j < nquad2; ++j)
1967  {
1968  for(i = 0; i < nquad1; ++i)
1969  {
1970  Blas::Dscal(nquad0,w1[i], &outarray[0]+i*nquad0 +
1971  j*nquad0*nquad1,1);
1972  }
1973  }
1974 
1975  // Multiply by integration constants in z-direction; need to
1976  // incorporate factor [(1-eta_3)/2]^2 into weights, but only if
1977  // using GLL quadrature points.
1978  switch(m_base[2]->GetPointsType())
1979  {
1980  // (2,0) Jacobi inner product.
1982  for(i = 0; i < nquad2; ++i)
1983  {
1984  Blas::Dscal(nquad0*nquad1, 0.25*w2[i],
1985  &outarray[0]+i*nquad0*nquad1, 1);
1986  }
1987  break;
1988 
1989  default:
1990  for(i = 0; i < nquad2; ++i)
1991  {
1992  Blas::Dscal(nquad0*nquad1,0.125*(1-z2[i])*(1-z2[i])*w2[i],
1993  &outarray[0]+i*nquad0*nquad1,1);
1994  }
1995  break;
1996  }
1997  }
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:58
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:216
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:169
int Nektar::StdRegions::StdPyrExp::v_NumBndryCoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 896 of file StdPyrExp.cpp.

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

Referenced by v_GetBoundaryMap(), and v_GetInteriorMap().

897  {
900  "BasisType is not a boundary interior form");
903  "BasisType is not a boundary interior form");
906  "BasisType is not a boundary interior form");
907 
908  int P = m_base[0]->GetNumModes();
909  int Q = m_base[1]->GetNumModes();
910  int R = m_base[2]->GetNumModes();
911 
914  }
Principle Modified Functions .
Definition: BasisType.h:51
Principle Modified Functions .
Definition: BasisType.h:49
int getNumberOfBndCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:266
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
Lagrange for SEM basis .
Definition: BasisType.h:53
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
Array< OneD, LibUtilities::BasisSharedPtr > m_base
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 253 of file StdPyrExp.cpp.

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

Referenced by v_PhysDeriv(), and v_StdPhysDeriv().

258  {
259  // PhysDerivative implementation based on Spen's book page 152.
260  int Qx = m_base[0]->GetNumPoints();
261  int Qy = m_base[1]->GetNumPoints();
262  int Qz = m_base[2]->GetNumPoints();
263 
264  Array<OneD, NekDouble> dEta_bar1(Qx*Qy*Qz,0.0);
265  Array<OneD, NekDouble> dXi2 (Qx*Qy*Qz,0.0);
266  Array<OneD, NekDouble> dEta3 (Qx*Qy*Qz,0.0);
267  PhysTensorDeriv(u_physical, dEta_bar1, dXi2, dEta3);
268 
269  Array<OneD, const NekDouble> eta_x, eta_y, eta_z;
270  eta_x = m_base[0]->GetZ();
271  eta_y = m_base[1]->GetZ();
272  eta_z = m_base[2]->GetZ();
273 
274  int i, j, k, n;
275 
276  for (k = 0, n = 0; k < Qz; ++k)
277  {
278  for (j = 0; j < Qy; ++j)
279  {
280  for (i = 0; i < Qx; ++i, ++n)
281  {
282  if (out_dxi1.num_elements() > 0)
283  out_dxi1[n] = 2.0/(1.0 - eta_z[k]) * dEta_bar1[n];
284  if (out_dxi2.num_elements() > 0)
285  out_dxi2[n] = 2.0/(1.0 - eta_z[k]) * dXi2[n];
286  if (out_dxi3.num_elements() > 0)
287  out_dxi3[n] = (1.0+eta_x[i])/(1.0-eta_z[k])*dEta_bar1[n] +
288  (1.0+eta_y[j])/(1.0-eta_z[k])*dXi2[n] + dEta3[n];
289  }
290  }
291  }
292  }
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...
Array< OneD, LibUtilities::BasisSharedPtr > m_base
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 294 of file StdPyrExp.cpp.

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

297  {
298  switch(dir)
299  {
300  case 0:
301  {
302  v_PhysDeriv(inarray, outarray, NullNekDouble1DArray,
304  break;
305  }
306 
307  case 1:
308  {
309  v_PhysDeriv(inarray, NullNekDouble1DArray, outarray,
311  break;
312  }
313 
314  case 2:
315  {
317  NullNekDouble1DArray, outarray);
318  break;
319  }
320 
321  default:
322  {
323  ASSERTL1(false,"input dir is out of range");
324  }
325  break;
326  }
327  }
static Array< OneD, NekDouble > NullNekDouble1DArray
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
Calculate the derivative of the physical points.
Definition: StdPyrExp.cpp:253
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
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 329 of file StdPyrExp.cpp.

References v_PhysDeriv().

333  {
334  StdPyrExp::v_PhysDeriv(inarray, out_d0, out_d1, out_d2);
335  }
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
Calculate the derivative of the physical points.
Definition: StdPyrExp.cpp:253
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 337 of file StdPyrExp.cpp.

References v_PhysDeriv().

340  {
341  StdPyrExp::v_PhysDeriv(dir, inarray, outarray);
342  }
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
Calculate the derivative of the physical points.
Definition: StdPyrExp.cpp:253

Member Data Documentation

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

Definition at line 256 of file StdPyrExp.h.

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

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

Definition at line 255 of file StdPyrExp.h.

Referenced by StdPyrExp().