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

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

#include <StdPrismExp.h>

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

Public Member Functions

 StdPrismExp ()
 
 StdPrismExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
 
 StdPrismExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc, NekDouble *coeffs, NekDouble *phys)
 
 StdPrismExp (const StdPrismExp &T)
 
 ~StdPrismExp ()
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion3D
 StdExpansion3D ()
 
 StdExpansion3D (int numcoeffs, const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
 
 StdExpansion3D (const StdExpansion3D &T)
 
virtual ~StdExpansion3D ()
 
void PhysTensorDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d1, Array< OneD, NekDouble > &outarray_d2, Array< OneD, NekDouble > &outarray_d3)
 Calculate the 3D derivative in the local tensor/collapsed coordinate at the physical points. More...
 
void BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
 
void IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
 
- 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
 
boost::shared_ptr< StdExpansionGetLinStdExp (void) const
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion ()
 
bool IsNodalNonTensorialExp ()
 
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space. More...
 
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space. More...
 
void 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)
 
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 GetFaceNumModes (const int fid, const Orientation faceOrient, int &numModes0, int &numModes1)
 
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_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)
 
virtual void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
 
virtual void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Forward transform from physical quadrature space stored in inarray and evaluate the expansion coefficients and store in outarray. More...
 
virtual void v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Calculate the inner product of inarray with respect to the basis B=base0*base1*base2 and put into outarray: More...
 
virtual void v_IProductWRTBase_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
virtual void v_IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
 
virtual void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Inner product of inarray over region with respect to the object's default expansion basis; output in outarray. More...
 
virtual void v_IProductWRTDerivBase_MatOp (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 
virtual void v_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 void v_GetFaceNumModes (const int fid, const Orientation faceOrient, int &numModes0, int &numModes1)
 
virtual int v_GetNverts () const
 
virtual int v_GetNedges () const
 
virtual int v_GetNfaces () const
 
virtual LibUtilities::ShapeType v_DetShapeType () const
 Return Shape of region, using ShapeType enum list; i.e. prism. More...
 
virtual int v_NumBndryCoeffs () const
 
virtual int v_NumDGBndryCoeffs () const
 
virtual int v_GetEdgeNcoeffs (const int i) const
 
virtual int v_GetTotalEdgeIntNcoeffs () const
 
virtual int v_GetFaceNcoeffs (const int i) const
 
virtual int v_GetFaceNumPoints (const int i) const
 
virtual LibUtilities::PointsKey v_GetFacePointsKey (const int i, const int j) const
 
virtual const
LibUtilities::BasisKey 
v_DetFaceBasisKey (const int i, const int k) const
 
virtual int v_GetFaceIntNcoeffs (const int i) const
 
virtual int v_GetTotalFaceIntNcoeffs () const
 
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
virtual LibUtilities::BasisType v_GetEdgeBasisType (const int i) const
 
virtual bool v_IsBoundaryInteriorExpansion ()
 
virtual void v_GetFaceToElementMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int P=-1, int Q=-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)
 
virtual void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void v_SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
virtual void v_ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion3D
virtual NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
virtual NekDouble v_PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 
virtual void v_LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
 
virtual void v_HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
 
virtual NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray)
 Integrates the specified function over the domain. More...
 
virtual void v_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 Member Functions

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

Additional Inherited Members

- Protected Attributes inherited from Nektar::StdRegions::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

Class representing a prismatic element in reference space.

Definition at line 49 of file StdPrismExp.h.

Constructor & Destructor Documentation

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

Definition at line 46 of file StdPrismExp.cpp.

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

Definition at line 50 of file StdPrismExp.cpp.

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

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  ASSERTL0(Ba.GetNumModes() <= Bc.GetNumModes(),
65  "order in 'a' direction is higher than order in 'c' direction");
66  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
int getNumberOfCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:286
StdExpansion()
Default Constructor.
Nektar::StdRegions::StdPrismExp::StdPrismExp ( const LibUtilities::BasisKey Ba,
const LibUtilities::BasisKey Bb,
const LibUtilities::BasisKey Bc,
NekDouble coeffs,
NekDouble phys 
)
Nektar::StdRegions::StdPrismExp::StdPrismExp ( const StdPrismExp T)

Definition at line 68 of file StdPrismExp.cpp.

69  : StdExpansion(T),
71  {
72  }
StdExpansion()
Default Constructor.
Nektar::StdRegions::StdPrismExp::~StdPrismExp ( )

Definition at line 76 of file StdPrismExp.cpp.

77  {
78  }

Member Function Documentation

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

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

Modes are numbered with the r index travelling fastest, followed by q and then p, and each q-r plane is of size (R+1-p). For example, with P=1, Q=2, R=3, the indexing inside each q-r plane (with r increasing upwards and q to the right) is:

p = 0: p = 1:

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

Note that in this element, we must have that $ P <= R $.

Definition at line 1932 of file StdPrismExp.cpp.

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

Referenced by v_GetBoundaryMap(), v_GetEdgeInteriorMap(), v_GetFaceInteriorMap(), v_GetFaceToElementMap(), v_GetInteriorMap(), v_GetVertexMap(), and v_IProductWRTBase_SumFacKernel().

1933  {
1934  int Q = m_base[1]->GetNumModes() - 1;
1935  int R = m_base[2]->GetNumModes() - 1;
1936 
1937  return r + // Skip along stacks (r-direction)
1938  q*(R+1-p) + // Skip along columns (q-direction)
1939  (Q+1)*(p*R + 1-(p-2)*(p-1)/2); // Skip along rows (p-direction)
1940  }
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdPrismExp::v_BwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual
Note
'r' (base[2]) runs fastest in this element.

Perform backwards transformation at the quadrature points:

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

In the prism this expansion becomes:

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

And sumfactorizing step of the form is as:\

$ f_{pr} (\xi_{3k}) = \sum_{r=0}^{Q_z} \hat u_{pqr} \psi_{pr}^b (\xi_{3k}),\\ g_{p} (\xi_{2j}, \xi_{3k}) = \sum_{r=0}^{Q_y} \psi_{p}^a (\xi_{2j}) f_{pr} (\xi_{3k}),\ \ u(\xi_{1i}, \xi_{2j}, \xi_{3k}) = \sum_{p=0}^{Q_x} \psi_{p}^a (\xi_{1i}) g_{p} (\xi_{2j}, \xi_{3k}). $

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalPrismExp.

Definition at line 242 of file StdPrismExp.cpp.

References ASSERTL1, Nektar::LibUtilities::eModified_B, Nektar::LibUtilities::eModified_C, Nektar::LibUtilities::eOrtho_B, Nektar::LibUtilities::eOrtho_C, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::GetNumPoints(), Nektar::StdRegions::StdExpansion::m_base, v_BwdTrans_SumFac(), and Vmath::Vcopy().

Referenced by v_FillMode().

244  {
247  "Basis[1] is not a general tensor type");
248 
251  "Basis[2] is not a general tensor type");
252 
253  if(m_base[0]->Collocation() &&
254  m_base[1]->Collocation() &&
255  m_base[2]->Collocation())
256  {
258  m_base[1]->GetNumPoints()*
259  m_base[2]->GetNumPoints(),
260  inarray, 1, outarray, 1);
261  }
262  else
263  {
264  StdPrismExp::v_BwdTrans_SumFac(inarray,outarray);
265  }
266  }
Principle Modified Functions .
Definition: BasisType.h:51
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:229
Principle Orthogonal Functions .
Definition: BasisType.h:47
Principle Modified Functions .
Definition: BasisType.h:50
Principle Orthogonal Functions .
Definition: BasisType.h:48
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
void Nektar::StdRegions::StdPrismExp::v_BwdTrans_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalPrismExp.

Definition at line 268 of file StdPrismExp.cpp.

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

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

270  {
271  int nquad1 = m_base[1]->GetNumPoints();
272  int nquad2 = m_base[2]->GetNumPoints();
273  int order0 = m_base[0]->GetNumModes();
274  int order1 = m_base[1]->GetNumModes();
275 
276  Array<OneD, NekDouble> wsp(nquad2*order1*order0 +
277  nquad1*nquad2*order0);
278 
279  BwdTrans_SumFacKernel(m_base[0]->GetBdata(),
280  m_base[1]->GetBdata(),
281  m_base[2]->GetBdata(),
282  inarray,outarray,wsp,true,true,true);
283  }
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)
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdPrismExp::v_BwdTrans_SumFacKernel ( const Array< OneD, const NekDouble > &  base0,
const Array< OneD, const NekDouble > &  base1,
const Array< OneD, const NekDouble > &  base2,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
Array< OneD, NekDouble > &  wsp,
bool  doCheckCollDir0,
bool  doCheckCollDir1,
bool  doCheckCollDir2 
)
protectedvirtual

Implements Nektar::StdRegions::StdExpansion3D.

Definition at line 286 of file StdPrismExp.cpp.

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

296  {
297  int i, mode;
298  int nquad0 = m_base[0]->GetNumPoints();
299  int nquad1 = m_base[1]->GetNumPoints();
300  int nquad2 = m_base[2]->GetNumPoints();
301  int nummodes0 = m_base[0]->GetNumModes();
302  int nummodes1 = m_base[1]->GetNumModes();
303  int nummodes2 = m_base[2]->GetNumModes();
304  Array<OneD, NekDouble> tmp0 = wsp;
305  Array<OneD, NekDouble> tmp1 = tmp0 + nquad2*nummodes1*nummodes0;
306 
307  for (i = mode = 0; i < nummodes0; ++i)
308  {
309  Blas::Dgemm('N', 'N', nquad2, nummodes1, nummodes2-i,
310  1.0, base2.get() + mode*nquad2, nquad2,
311  inarray.get() + mode*nummodes1, nummodes2-i,
312  0.0, tmp0.get() + i*nquad2*nummodes1, nquad2);
313  mode += nummodes2-i;
314  }
315 
317  {
318  for(i = 0; i < nummodes1; i++)
319  {
320  Blas::Daxpy(nquad2,inarray[1+i*nummodes2],base2.get()+nquad2,1,
321  tmp0.get()+nquad2*(nummodes1+i),1);
322  }
323  }
324 
325  for (i = 0; i < nummodes0; i++)
326  {
327  Blas::Dgemm('N', 'T', nquad1, nquad2, nummodes1,
328  1.0, base1.get(), nquad1,
329  tmp0.get() + i*nquad2*nummodes1, nquad2,
330  0.0, tmp1.get() + i*nquad2*nquad1, nquad1);
331  }
332 
333  Blas::Dgemm('N', 'T', nquad0, nquad2*nquad1, nummodes0,
334  1.0, base0.get(), nquad0,
335  tmp1.get(), nquad2*nquad1,
336  0.0, outarray.get(), nquad0);
337  }
Principle Modified Functions .
Definition: BasisType.h:49
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
Array< OneD, LibUtilities::BasisSharedPtr > m_base
int Nektar::StdRegions::StdPrismExp::v_CalcNumberOfCoefficients ( const std::vector< unsigned int > &  nummodes,
int &  modes_offset 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1000 of file StdPrismExp.cpp.

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

1002  {
1004  nummodes[modes_offset],
1005  nummodes[modes_offset+1],
1006  nummodes[modes_offset+2]);
1007 
1008  modes_offset += 3;
1009  return nmodes;
1010  }
int getNumberOfCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:286
DNekMatSharedPtr Nektar::StdRegions::StdPrismExp::v_CreateStdMatrix ( const StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 1907 of file StdPrismExp.cpp.

References v_GenMatrix().

1908  {
1909  return v_GenMatrix(mkey);
1910  }
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey)
const LibUtilities::BasisKey Nektar::StdRegions::StdPrismExp::v_DetFaceBasisKey ( const int  i,
const int  k 
) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 961 of file StdPrismExp.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().

963  {
964  ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
965  ASSERTL2(k >= 0 && k <= 1, "basis key id is out of range");
966 
967  switch(i)
968  {
969  case 0:
970  {
971  return EvaluateQuadFaceBasisKey(k,
972  m_base[k]->GetBasisType(),
973  m_base[k]->GetNumPoints(),
974  m_base[k]->GetNumModes());
975  }
976  case 2:
977  case 4:
978  {
979  return EvaluateQuadFaceBasisKey(k,
980  m_base[k+1]->GetBasisType(),
981  m_base[k+1]->GetNumPoints(),
982  m_base[k+1]->GetNumModes());
983  }
984  case 1:
985  case 3:
986  {
987  return EvaluateTriFaceBasisKey(k,
988  m_base[2*k]->GetBasisType(),
989  m_base[2*k]->GetNumPoints(),
990  m_base[2*k]->GetNumModes());
991 
992  }
993  break;
994  }
995 
996  // Should never get here.
998  }
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:250
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::StdPrismExp::v_DetShapeType ( ) const
protectedvirtual

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 796 of file StdPrismExp.cpp.

References Nektar::LibUtilities::ePrism.

797  {
798  return LibUtilities::ePrism;
799  }
void Nektar::StdRegions::StdPrismExp::v_FillMode ( const int  mode,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalPrismExp.

Definition at line 724 of file StdPrismExp.cpp.

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

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

725  {
726  Array<OneD, NekDouble> tmp(m_ncoeffs,0.0);
727  tmp[mode] = 1.0;
728  StdPrismExp::v_BwdTrans(tmp, outarray);
729  }
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void Nektar::StdRegions::StdPrismExp::v_FwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

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

Inputs:

  • inarray: array of physical quadrature points to be transformed

Outputs:

  • outarray: updated array of expansion coefficients.

Implements Nektar::StdRegions::StdExpansion.

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

Definition at line 350 of file StdPrismExp.cpp.

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

352  {
353  v_IProductWRTBase(inarray, outarray);
354 
355  // Get Mass matrix inverse
356  StdMatrixKey masskey(eInvMass,DetShapeType(),*this);
357  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
358 
359  // copy inarray in case inarray == outarray
360  DNekVec in (m_ncoeffs, outarray);
361  DNekVec out(m_ncoeffs, outarray, eWrapper);
362 
363  out = (*matsys)*in;
364  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:470
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Calculate the inner product of inarray with respect to the basis B=base0*base1*base2 and put into out...
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:706
NekVector< NekDouble > DNekVec
Definition: NekTypeDefs.hpp:49
DNekMatSharedPtr Nektar::StdRegions::StdPrismExp::v_GenMatrix ( const StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 1819 of file StdPrismExp.cpp.

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

Referenced by v_CreateStdMatrix().

1820  {
1821 
1822  MatrixType mtype = mkey.GetMatrixType();
1823 
1824  DNekMatSharedPtr Mat;
1825 
1826  switch(mtype)
1827  {
1829  {
1830  int nq0 = m_base[0]->GetNumPoints();
1831  int nq1 = m_base[1]->GetNumPoints();
1832  int nq2 = m_base[2]->GetNumPoints();
1833  int nq;
1834 
1835  // take definition from key
1836  if(mkey.ConstFactorExists(eFactorConst))
1837  {
1838  nq = (int) mkey.GetConstFactor(eFactorConst);
1839  }
1840  else
1841  {
1842  nq = max(nq0,max(nq1,nq2));
1843  }
1844 
1846  getNumberOfCoefficients (nq,nq,nq);
1847  Array<OneD, Array<OneD, NekDouble> > coords (neq);
1848  Array<OneD, NekDouble> coll (3);
1849  Array<OneD, DNekMatSharedPtr> I (3);
1850  Array<OneD, NekDouble> tmp (nq0);
1851 
1853  AllocateSharedPtr(neq,nq0*nq1*nq2);
1854  int cnt = 0;
1855  for(int i = 0; i < nq; ++i)
1856  {
1857  for(int j = 0; j < nq; ++j)
1858  {
1859  for(int k = 0; k < nq-i; ++k,++cnt)
1860  {
1861  coords[cnt] = Array<OneD, NekDouble>(3);
1862  coords[cnt][0] = -1.0 + 2*k/(NekDouble)(nq-1);
1863  coords[cnt][1] = -1.0 + 2*j/(NekDouble)(nq-1);
1864  coords[cnt][2] = -1.0 + 2*i/(NekDouble)(nq-1);
1865  }
1866  }
1867  }
1868 
1869  for(int i = 0; i < neq; ++i)
1870  {
1871  LocCoordToLocCollapsed(coords[i],coll);
1872 
1873  I[0] = m_base[0]->GetI(coll );
1874  I[1] = m_base[1]->GetI(coll+1);
1875  I[2] = m_base[2]->GetI(coll+2);
1876 
1877  // interpolate first coordinate direction
1878  NekDouble fac;
1879  for( int k = 0; k < nq2; ++k)
1880  {
1881  for (int j = 0; j < nq1; ++j)
1882  {
1883 
1884  fac = (I[1]->GetPtr())[j]*(I[2]->GetPtr())[k];
1885  Vmath::Smul(nq0,fac,I[0]->GetPtr(),1,tmp,1);
1886 
1887  Vmath::Vcopy(nq0, &tmp[0], 1,
1888  Mat->GetRawPtr() +
1889  k * nq0 * nq1 * neq +
1890  j * nq0 * neq + i,
1891  neq);
1892  }
1893  }
1894  }
1895  }
1896  break;
1897  default:
1898  {
1900  }
1901  break;
1902  }
1903 
1904  return Mat;
1905  }
int getNumberOfCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:286
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
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:213
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
double NekDouble
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
void Nektar::StdRegions::StdPrismExp::v_GetBoundaryMap ( Array< OneD, unsigned int > &  outarray)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalPrismExp.

Definition at line 1753 of file StdPrismExp.cpp.

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

1754  {
1757  "BasisType is not a boundary interior form");
1760  "BasisType is not a boundary interior form");
1763  "BasisType is not a boundary interior form");
1764 
1765  int P = m_base[0]->GetNumModes() - 1, p;
1766  int Q = m_base[1]->GetNumModes() - 1, q;
1767  int R = m_base[2]->GetNumModes() - 1, r;
1768  int idx = 0;
1769 
1770  int nBnd = NumBndryCoeffs();
1771 
1772  if (maparray.num_elements() != nBnd)
1773  {
1774  maparray = Array<OneD, unsigned int>(nBnd);
1775  }
1776 
1777  // Loop over all boundary modes (in ascending order).
1778  for (p = 0; p <= P; ++p)
1779  {
1780  // First two q-r planes are entirely boundary modes.
1781  if (p <= 1)
1782  {
1783  for (q = 0; q <= Q; ++q)
1784  {
1785  for (r = 0; r <= R-p; ++r)
1786  {
1787  maparray[idx++] = GetMode(p,q,r);
1788  }
1789  }
1790  }
1791  else
1792  {
1793  // Remaining q-r planes contain boundary modes on the two
1794  // left-hand sides and bottom edge.
1795  for (q = 0; q <= Q; ++q)
1796  {
1797  if (q <= 1)
1798  {
1799  for (r = 0; r <= R-p; ++r)
1800  {
1801  maparray[idx++] = GetMode(p,q,r);
1802  }
1803  }
1804  else
1805  {
1806  maparray[idx++] = GetMode(p,q,0);
1807  }
1808  }
1809  }
1810  }
1811  }
Principle Modified Functions .
Definition: BasisType.h:49
Principle Modified Functions .
Definition: BasisType.h:50
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
int GetMode(int I, int J, int K)
Compute the local mode number in the expansion for a particular tensorial combination.
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:228
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdPrismExp::v_GetCoords ( Array< OneD, NekDouble > &  xi_x,
Array< OneD, NekDouble > &  xi_y,
Array< OneD, NekDouble > &  xi_z 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PrismExp.

Definition at line 700 of file StdPrismExp.cpp.

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

703  {
704  Array<OneD, const NekDouble> etaBar_x = m_base[0]->GetZ();
705  Array<OneD, const NekDouble> eta_y = m_base[1]->GetZ();
706  Array<OneD, const NekDouble> eta_z = m_base[2]->GetZ();
707  int Qx = GetNumPoints(0);
708  int Qy = GetNumPoints(1);
709  int Qz = GetNumPoints(2);
710 
711  // Convert collapsed coordinates into cartesian coordinates: eta --> xi
712  for (int k = 0; k < Qz; ++k) {
713  for (int j = 0; j < Qy; ++j) {
714  for (int i = 0; i < Qx; ++i) {
715  int s = i + Qx*(j + Qy*k);
716  xi_x[s] = (1.0 - eta_z[k])*(1.0 + etaBar_x[i]) / 2.0 - 1.0;
717  xi_y[s] = eta_y[j];
718  xi_z[s] = eta_z[k];
719  }
720  }
721  }
722  }
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::StdPrismExp::v_GetEdgeBasisType ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1012 of file StdPrismExp.cpp.

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

1013  {
1014  ASSERTL2(i >= 0 && i <= 8, "edge id is out of range");
1015  if (i == 0 || i == 2)
1016  {
1017  return GetBasisType(0);
1018  }
1019  else if (i == 1 || i == 3 || i == 8)
1020  {
1021  return GetBasisType(1);
1022  }
1023  else
1024  {
1025  return GetBasisType(2);
1026  }
1027  }
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:250
void Nektar::StdRegions::StdPrismExp::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 1419 of file StdPrismExp.cpp.

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

1424  {
1425  int i;
1426  bool signChange;
1427  const int P = m_base[0]->GetNumModes() - 1;
1428  const int Q = m_base[1]->GetNumModes() - 1;
1429  const int R = m_base[2]->GetNumModes() - 1;
1430  const int nEdgeIntCoeffs = v_GetEdgeNcoeffs(eid) - 2;
1431 
1432  if (maparray.num_elements() != nEdgeIntCoeffs)
1433  {
1434  maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
1435  }
1436 
1437  if(signarray.num_elements() != nEdgeIntCoeffs)
1438  {
1439  signarray = Array<OneD, int>(nEdgeIntCoeffs,1);
1440  }
1441  else
1442  {
1443  fill(signarray.get(), signarray.get()+nEdgeIntCoeffs, 1);
1444  }
1445 
1446  // If edge is oriented backwards, change sign of modes which have
1447  // degree 2n+1, n >= 1.
1448  signChange = edgeOrient == eBackwards;
1449 
1450  switch (eid)
1451  {
1452  case 0:
1453  for (i = 2; i <= P; ++i)
1454  {
1455  maparray[i-2] = GetMode(i,0,0);
1456  }
1457  break;
1458 
1459  case 1:
1460  for (i = 2; i <= Q; ++i)
1461  {
1462  maparray[i-2] = GetMode(1,i,0);
1463  }
1464  break;
1465 
1466  case 2:
1467  // Base quad; reverse direction.
1468  //signChange = !signChange;
1469  for (i = 2; i <= P; ++i)
1470  {
1471  maparray[i-2] = GetMode(i,1,0);
1472  }
1473  break;
1474 
1475  case 3:
1476  // Base quad; reverse direction.
1477  //signChange = !signChange;
1478  for (i = 2; i <= Q; ++i)
1479  {
1480  maparray[i-2] = GetMode(0,i,0);
1481  }
1482  break;
1483 
1484  case 4:
1485  for (i = 2; i <= R; ++i)
1486  {
1487  maparray[i-2] = GetMode(0,0,i);
1488  }
1489  break;
1490 
1491  case 5:
1492  for (i = 1; i <= R-1; ++i)
1493  {
1494  maparray[i-1] = GetMode(1,0,i);
1495  }
1496  break;
1497 
1498  case 6:
1499  for (i = 1; i <= R-1; ++i)
1500  {
1501  maparray[i-1] = GetMode(1,1,i);
1502  }
1503  break;
1504 
1505  case 7:
1506  for (i = 2; i <= R; ++i)
1507  {
1508  maparray[i-2] = GetMode(0,1,i);
1509  }
1510  break;
1511 
1512  case 8:
1513  for (i = 2; i <= Q; ++i)
1514  {
1515  maparray[i-2] = GetMode(0,i,1);
1516  }
1517  break;
1518 
1519  default:
1520  ASSERTL0(false, "Edge not defined.");
1521  break;
1522  }
1523 
1524  if (signChange)
1525  {
1526  for (i = 1; i < nEdgeIntCoeffs; i += 2)
1527  {
1528  signarray[i] = -1;
1529  }
1530  }
1531  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
virtual int v_GetEdgeNcoeffs(const int i) const
int GetMode(int I, int J, int K)
Compute the local mode number in the expansion for a particular tensorial combination.
Array< OneD, LibUtilities::BasisSharedPtr > m_base
int Nektar::StdRegions::StdPrismExp::v_GetEdgeNcoeffs ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 842 of file StdPrismExp.cpp.

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

Referenced by v_GetEdgeInteriorMap().

843  {
844  ASSERTL2(i >= 0 && i <= 8, "edge id is out of range");
845 
846  if (i == 0 || i == 2)
847  {
848  return GetBasisNumModes(0);
849  }
850  else if (i == 1 || i == 3 || i == 8)
851  {
852  return GetBasisNumModes(1);
853  }
854  else
855  {
856  return GetBasisNumModes(2);
857  }
858  }
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:250
void Nektar::StdRegions::StdPrismExp::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 1533 of file StdPrismExp.cpp.

References ASSERTL0, GetMode(), Nektar::StdRegions::StdExpansion::m_base, class_topology::P, CellMLToNektar.cellml_metadata::p, and v_GetFaceIntNcoeffs().

1538  {
1539  const int P = m_base[0]->GetNumModes() - 1;
1540  const int Q = m_base[1]->GetNumModes() - 1;
1541  const int R = m_base[2]->GetNumModes() - 1;
1542  const int nFaceIntCoeffs = v_GetFaceIntNcoeffs(fid);
1543  int p, q, r, idx = 0;
1544  int nummodesA = 0;
1545  int nummodesB = 0;
1546  int i = 0;
1547  int j = 0;
1548 
1549  if (maparray.num_elements() != nFaceIntCoeffs)
1550  {
1551  maparray = Array<OneD, unsigned int>(nFaceIntCoeffs);
1552  }
1553 
1554  if (signarray.num_elements() != nFaceIntCoeffs)
1555  {
1556  signarray = Array<OneD, int>(nFaceIntCoeffs, 1);
1557  }
1558  else
1559  {
1560  fill(signarray.get(), signarray.get()+nFaceIntCoeffs, 1);
1561  }
1562 
1563  // Set up an array indexing for quad faces, since the ordering may
1564  // need to be transposed depending on orientation.
1565  Array<OneD, int> arrayindx(nFaceIntCoeffs);
1566  if (fid != 1 && fid != 3)
1567  {
1568  if (fid == 0) // Base quad
1569  {
1570  nummodesA = P-1;
1571  nummodesB = Q-1;
1572  }
1573  else // front and back quad
1574  {
1575  nummodesA = Q-1;
1576  nummodesB = R-1;
1577  }
1578 
1579  for (i = 0; i < nummodesB; i++)
1580  {
1581  for (j = 0; j < nummodesA; j++)
1582  {
1583  if (faceOrient < 9)
1584  {
1585  arrayindx[i*nummodesA+j] = i*nummodesA+j;
1586  }
1587  else
1588  {
1589  arrayindx[i*nummodesA+j] = j*nummodesB+i;
1590  }
1591  }
1592  }
1593  }
1594 
1595  switch (fid)
1596  {
1597  case 0: // Bottom quad
1598  for (q = 2; q <= Q; ++q)
1599  {
1600  for (p = 2; p <= P; ++p)
1601  {
1602  maparray[arrayindx[(q-2)*nummodesA+(p-2)]] = GetMode(p,q,0);
1603  }
1604  }
1605  break;
1606 
1607  case 1: // Left triangle
1608  for (p = 2; p <= P; ++p)
1609  {
1610  for (r = 1; r <= R-p; ++r)
1611  {
1612  if ((int)faceOrient == 7)
1613  {
1614  signarray[idx] = p % 2 ? -1 : 1;
1615  }
1616  maparray[idx++] = GetMode(p,0,r);
1617  }
1618  }
1619  break;
1620 
1621  case 2: // Slanted quad
1622  for (r = 1; r <= R-1; ++r)
1623  {
1624  for (q = 2; q <= Q; ++q)
1625  {
1626  maparray[arrayindx[(r-1)*nummodesA+(q-2)]] = GetMode(1, q, r);
1627  }
1628  }
1629  break;
1630 
1631  case 3: // Right triangle
1632  for (p = 2; p <= P; ++p)
1633  {
1634  for (r = 1; r <= R-p; ++r)
1635  {
1636  if ((int)faceOrient == 7)
1637  {
1638  signarray[idx] = p % 2 ? -1 : 1;
1639  }
1640  maparray[idx++] = GetMode(p, 1, r);
1641  }
1642  }
1643  break;
1644 
1645  case 4: // Back quad
1646  for (r = 2; r <= R; ++r)
1647  {
1648  for (q = 2; q <= Q; ++q)
1649  {
1650  maparray[arrayindx[(r-2)*nummodesA+(q-2)]] = GetMode(0, q, r);
1651  }
1652  }
1653  break;
1654 
1655  default:
1656  ASSERTL0(false, "Face interior map not available.");
1657  }
1658 
1659  // Triangular faces are processed in the above switch loop; for
1660  // remaining quad faces, set up orientation if necessary.
1661  if (fid == 1 || fid == 3)
1662  return;
1663 
1664  if (faceOrient == 6 || faceOrient == 8 ||
1665  faceOrient == 11 || faceOrient == 12)
1666  {
1667  if (faceOrient < 9)
1668  {
1669  for (i = 1; i < nummodesB; i += 2)
1670  {
1671  for (j = 0; j < nummodesA; j++)
1672  {
1673  signarray[arrayindx[i*nummodesA+j]] *= -1;
1674  }
1675  }
1676  }
1677  else
1678  {
1679  for (i = 0; i < nummodesB; i++)
1680  {
1681  for (j = 1; j < nummodesA; j += 2)
1682  {
1683  signarray[arrayindx[i*nummodesA+j]] *= -1;
1684  }
1685  }
1686  }
1687  }
1688 
1689  if (faceOrient == 7 || faceOrient == 8 ||
1690  faceOrient == 10 || faceOrient == 12)
1691  {
1692  if (faceOrient < 9)
1693  {
1694  for (i = 0; i < nummodesB; i++)
1695  {
1696  for (j = 1; j < nummodesA; j += 2)
1697  {
1698  signarray[arrayindx[i*nummodesA+j]] *= -1;
1699  }
1700  }
1701  }
1702  else
1703  {
1704  for (i = 1; i < nummodesB; i += 2)
1705  {
1706  for (j = 0; j < nummodesA; j++)
1707  {
1708  signarray[arrayindx[i*nummodesA+j]] *= -1;
1709  }
1710  }
1711  }
1712  }
1713  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
virtual int v_GetFaceIntNcoeffs(const int i) const
int GetMode(int I, int J, int K)
Compute the local mode number in the expansion for a particular tensorial combination.
Array< OneD, LibUtilities::BasisSharedPtr > m_base
int Nektar::StdRegions::StdPrismExp::v_GetFaceIntNcoeffs ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 887 of file StdPrismExp.cpp.

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

Referenced by v_GetFaceInteriorMap().

888  {
889  ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
890 
891  int Pi = GetBasisNumModes(0) - 2;
892  int Qi = GetBasisNumModes(1) - 2;
893  int Ri = GetBasisNumModes(2) - 2;
894 
895  if (i == 0)
896  {
897  return Pi * Qi;
898  }
899  else if (i == 1 || i == 3)
900  {
901  return Pi * (2*Ri - Pi - 1) / 2;
902  }
903  else
904  {
905  return Qi * Ri;
906  }
907  }
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:250
int Nektar::StdRegions::StdPrismExp::v_GetFaceNcoeffs ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 869 of file StdPrismExp.cpp.

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

870  {
871  ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
872  if (i == 0)
873  {
874  return GetBasisNumModes(0)*GetBasisNumModes(1);
875  }
876  else if (i == 1 || i == 3)
877  {
878  int P = GetBasisNumModes(0)-1, Q = GetBasisNumModes(2)-1;
879  return Q+1 + (P*(1 + 2*Q - P))/2;
880  }
881  else
882  {
883  return GetBasisNumModes(1)*GetBasisNumModes(2);
884  }
885  }
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:250
void Nektar::StdRegions::StdPrismExp::v_GetFaceNumModes ( const int  fid,
const Orientation  faceOrient,
int &  numModes0,
int &  numModes1 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 731 of file StdPrismExp.cpp.

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

736  {
737  int nummodes [3] = {m_base[0]->GetNumModes(),
738  m_base[1]->GetNumModes(),
739  m_base[2]->GetNumModes()};
740  switch(fid)
741  {
742  // base quad
743  case 0:
744  {
745  numModes0 = nummodes[0];
746  numModes1 = nummodes[1];
747  }
748  break;
749  // front and back quad
750  case 2:
751  case 4:
752  {
753  numModes0 = nummodes[1];
754  numModes1 = nummodes[2];
755  }
756  break;
757  // triangles
758  case 1:
759  case 3:
760  {
761  numModes0 = nummodes[0];
762  numModes1 = nummodes[2];
763  }
764  break;
765  }
766 
767  if ( faceOrient >= 9 )
768  {
769  std::swap(numModes0, numModes1);
770  }
771  }
Array< OneD, LibUtilities::BasisSharedPtr > m_base
int Nektar::StdRegions::StdPrismExp::v_GetFaceNumPoints ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 920 of file StdPrismExp.cpp.

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

921  {
922  ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
923 
924  if (i == 0)
925  {
926  return m_base[0]->GetNumPoints()*
927  m_base[1]->GetNumPoints();
928  }
929  else if (i == 1 || i == 3)
930  {
931  return m_base[0]->GetNumPoints()*
932  m_base[2]->GetNumPoints();
933  }
934  else
935  {
936  return m_base[1]->GetNumPoints()*
937  m_base[2]->GetNumPoints();
938  }
939  }
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:250
Array< OneD, LibUtilities::BasisSharedPtr > m_base
LibUtilities::PointsKey Nektar::StdRegions::StdPrismExp::v_GetFacePointsKey ( const int  i,
const int  j 
) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 941 of file StdPrismExp.cpp.

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

943  {
944  ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
945  ASSERTL2(j == 0 || j == 1, "face direction is out of range");
946 
947  if (i == 0)
948  {
949  return m_base[j]->GetPointsKey();
950  }
951  else if (i == 1 || i == 3)
952  {
953  return m_base[2*j]->GetPointsKey();
954  }
955  else
956  {
957  return m_base[j+1]->GetPointsKey();
958  }
959  }
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:250
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdPrismExp::v_GetFaceToElementMap ( const int  fid,
const Orientation  faceOrient,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray,
int  P = -1,
int  Q = -1 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1041 of file StdPrismExp.cpp.

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

1048  {
1050  "Method only implemented if BasisType is identical"
1051  "in x and y directions");
1054  "Method only implemented for Modified_A BasisType"
1055  "(x and y direction) and Modified_B BasisType (z "
1056  "direction)");
1057 
1058  int i, j, k, p, q, r, nFaceCoeffs, idx = 0;
1059  int nummodesA=0, nummodesB=0;
1060 
1061  switch (fid)
1062  {
1063  case 0:
1064  nummodesA = m_base[0]->GetNumModes();
1065  nummodesB = m_base[1]->GetNumModes();
1066  break;
1067  case 1:
1068  case 3:
1069  nummodesA = m_base[0]->GetNumModes();
1070  nummodesB = m_base[2]->GetNumModes();
1071  break;
1072  case 2:
1073  case 4:
1074  nummodesA = m_base[1]->GetNumModes();
1075  nummodesB = m_base[2]->GetNumModes();
1076  break;
1077  }
1078 
1079  bool CheckForZeroedModes = false;
1080 
1081  if (P == -1)
1082  {
1083  P = nummodesA;
1084  Q = nummodesB;
1085  nFaceCoeffs = GetFaceNcoeffs(fid);
1086  }
1087  else if (fid == 1 || fid == 3)
1088  {
1089  nFaceCoeffs = P*(2*Q-P+1)/2;
1090  CheckForZeroedModes = true;
1091  }
1092  else
1093  {
1094  nFaceCoeffs = P*Q;
1095  CheckForZeroedModes = true;
1096  }
1097 
1098  // Allocate the map array and sign array; set sign array to ones (+)
1099  if (maparray.num_elements() != nFaceCoeffs)
1100  {
1101  maparray = Array<OneD, unsigned int>(nFaceCoeffs);
1102  }
1103 
1104  if (signarray.num_elements() != nFaceCoeffs)
1105  {
1106  signarray = Array<OneD, int>(nFaceCoeffs,1);
1107  }
1108  else
1109  {
1110  fill(signarray.get(), signarray.get() + nFaceCoeffs, 1);
1111  }
1112 
1113  // Set up an array indexing for quads, since the ordering may need
1114  // to be transposed.
1115  Array<OneD, int> arrayindx(nFaceCoeffs,-1);
1116 
1117  if (fid != 1 && fid != 3)
1118  {
1119  for (i = 0; i < Q; i++)
1120  {
1121  for (j = 0; j < P; j++)
1122  {
1123  if (faceOrient < 9)
1124  {
1125  arrayindx[i*P+j] = i*P+j;
1126  }
1127  else
1128  {
1129  arrayindx[i*P+j] = j*Q+i;
1130  }
1131  }
1132  }
1133  }
1134 
1135  // Set up ordering inside each 2D face. Also for triangular faces,
1136  // populate signarray.
1137  switch (fid)
1138  {
1139  case 0: // Bottom quad
1140  for (q = 0; q < Q; ++q)
1141  {
1142  for (p = 0; p < P; ++p)
1143  {
1144  maparray[arrayindx[q*P+p]] = GetMode(p,q,0);
1145  }
1146  }
1147  break;
1148 
1149  case 1: // Left triangle
1150  for (p = 0; p < P; ++p)
1151  {
1152  for (r = 0; r < Q-p; ++r)
1153  {
1154  if ((int)faceOrient == 7 && p > 1)
1155  {
1156  signarray[idx] = p % 2 ? -1 : 1;
1157  }
1158  maparray[idx++] = GetMode(p,0,r);
1159  }
1160  }
1161  break;
1162 
1163  case 2: // Slanted quad
1164  for (q = 0; q < P; ++q)
1165  {
1166  maparray[arrayindx[q]] = GetMode(1,q,0);
1167  }
1168  for (q = 0; q < P; ++q)
1169  {
1170  maparray[arrayindx[P+q]] = GetMode(0,q,1);
1171  }
1172  for (r = 1; r < Q-1; ++r)
1173  {
1174  for (q = 0; q < P; ++q)
1175  {
1176  maparray[arrayindx[(r+1)*P+q]] = GetMode(1,q,r);
1177  }
1178  }
1179  break;
1180 
1181  case 3: // Right triangle
1182  for (p = 0; p < P; ++p)
1183  {
1184  for (r = 0; r < Q-p; ++r)
1185  {
1186  if ((int)faceOrient == 7 && p > 1)
1187  {
1188  signarray[idx] = p % 2 ? -1 : 1;
1189  }
1190  maparray[idx++] = GetMode(p, 1, r);
1191  }
1192  }
1193  break;
1194 
1195  case 4: // Rear quad
1196  for (r = 0; r < Q; ++r)
1197  {
1198  for (q = 0; q < P; ++q)
1199  {
1200  maparray[arrayindx[r*P+q]] = GetMode(0, q, r);
1201  }
1202  }
1203  break;
1204 
1205  default:
1206  ASSERTL0(false, "Face to element map unavailable.");
1207  }
1208 
1209  if (fid == 1 || fid == 3)
1210  {
1211  if(CheckForZeroedModes)
1212  {
1213  // zero signmap and set maparray to zero if elemental
1214  // modes are not as large as face modesl
1215  idx = 0;
1216  for (j = 0; j < nummodesA; ++j)
1217  {
1218  idx += nummodesB-j;
1219  for (k = nummodesB-j; k < Q-j; ++k)
1220  {
1221  signarray[idx] = 0.0;
1222  maparray[idx++] = maparray[0];
1223  }
1224  }
1225 
1226  for (j = nummodesA; j < P; ++j)
1227  {
1228  for (k = 0; k < Q-j; ++k)
1229  {
1230  signarray[idx] = 0.0;
1231  maparray[idx++] = maparray[0];
1232  }
1233  }
1234  }
1235 
1236 
1237  // Triangles only have one possible orientation (base
1238  // direction reversed); swap edge modes.
1239  if ((int)faceOrient == 7)
1240  {
1241  swap(maparray[0], maparray[P]);
1242  for (i = 1; i < P-1; ++i)
1243  {
1244  swap(maparray[i+1], maparray[P+i]);
1245  }
1246  }
1247  }
1248  else
1249  {
1250 
1251  if(CheckForZeroedModes)
1252  {
1253  // zero signmap and set maparray to zero if elemental
1254  // modes are not as large as face modesl
1255  for (j = 0; j < nummodesA; ++j)
1256  {
1257  for (k = nummodesB; k < Q; ++k)
1258  {
1259  signarray[arrayindx[j+k*P]] = 0.0;
1260  maparray[arrayindx[j+k*P]] = maparray[0];
1261  }
1262  }
1263 
1264  for (j = nummodesA; j < P; ++j)
1265  {
1266  for (k = 0; k < Q; ++k)
1267  {
1268  signarray[arrayindx[j+k*P]] = 0.0;
1269  maparray[arrayindx[j+k*P]] = maparray[0];
1270  }
1271  }
1272  }
1273 
1274  // The code below is exactly the same as that taken from
1275  // StdHexExp and reverses the 'b' and 'a' directions as
1276  // appropriate (1st and 2nd if statements respectively) in
1277  // quadrilateral faces.
1278  if (faceOrient == 6 || faceOrient == 8 ||
1279  faceOrient == 11 || faceOrient == 12)
1280  {
1281  if (faceOrient < 9)
1282  {
1283  for (i = 3; i < Q; i += 2)
1284  {
1285  for (j = 0; j < P; j++)
1286  {
1287  signarray[arrayindx[i*P+j]] *= -1;
1288  }
1289  }
1290 
1291  for (i = 0; i < P; i++)
1292  {
1293  swap(maparray [i], maparray [i+P]);
1294  swap(signarray[i], signarray[i+P]);
1295  }
1296  }
1297  else
1298  {
1299  for (i = 0; i < Q; i++)
1300  {
1301  for (j = 3; j < P; j += 2)
1302  {
1303  signarray[arrayindx[i*P+j]] *= -1;
1304  }
1305  }
1306 
1307  for (i = 0; i < Q; i++)
1308  {
1309  swap (maparray [i], maparray [i+Q]);
1310  swap (signarray[i], signarray[i+Q]);
1311  }
1312  }
1313  }
1314 
1315  if (faceOrient == 7 || faceOrient == 8 ||
1316  faceOrient == 10 || faceOrient == 12)
1317  {
1318  if (faceOrient < 9)
1319  {
1320  for (i = 0; i < Q; i++)
1321  {
1322  for (j = 3; j < P; j += 2)
1323  {
1324  signarray[arrayindx[i*P+j]] *= -1;
1325  }
1326  }
1327 
1328  for(i = 0; i < Q; i++)
1329  {
1330  swap(maparray [i*P], maparray [i*P+1]);
1331  swap(signarray[i*P], signarray[i*P+1]);
1332  }
1333  }
1334  else
1335  {
1336  for (i = 3; i < Q; i += 2)
1337  {
1338  for (j = 0; j < P; j++)
1339  {
1340  signarray[arrayindx[i*P+j]] *= -1;
1341  }
1342  }
1343 
1344  for (i = 0; i < P; i++)
1345  {
1346  swap(maparray [i*Q], maparray [i*Q+1]);
1347  swap(signarray[i*Q], signarray[i*Q+1]);
1348  }
1349  }
1350  }
1351  }
1352  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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
Principle Modified Functions .
Definition: BasisType.h:50
int GetFaceNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th face.
Definition: StdExpansion.h:354
int GetMode(int I, int J, int K)
Compute the local mode number in the expansion for a particular tensorial combination.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdPrismExp::v_GetInteriorMap ( Array< OneD, unsigned int > &  outarray)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalPrismExp.

Definition at line 1715 of file StdPrismExp.cpp.

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

1716  {
1719  "BasisType is not a boundary interior form");
1722  "BasisType is not a boundary interior form");
1725  "BasisType is not a boundary interior form");
1726 
1727  int P = m_base[0]->GetNumModes() - 1, p;
1728  int Q = m_base[1]->GetNumModes() - 1, q;
1729  int R = m_base[2]->GetNumModes() - 1, r;
1730 
1731  int nIntCoeffs = m_ncoeffs - NumBndryCoeffs();
1732 
1733  if(outarray.num_elements()!=nIntCoeffs)
1734  {
1735  outarray = Array<OneD, unsigned int>(nIntCoeffs);
1736  }
1737 
1738  int idx = 0;
1739 
1740  // Loop over all interior modes.
1741  for (p = 2; p <= P; ++p)
1742  {
1743  for (q = 2; q <= Q; ++q)
1744  {
1745  for (r = 1; r <= R-p; ++r)
1746  {
1747  outarray[idx++] = GetMode(p,q,r);
1748  }
1749  }
1750  }
1751  }
Principle Modified Functions .
Definition: BasisType.h:49
Principle Modified Functions .
Definition: BasisType.h:50
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
int GetMode(int I, int J, int K)
Compute the local mode number in the expansion for a particular tensorial combination.
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:228
Array< OneD, LibUtilities::BasisSharedPtr > m_base
int Nektar::StdRegions::StdPrismExp::v_GetNedges ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 782 of file StdPrismExp.cpp.

783  {
784  return 9;
785  }
int Nektar::StdRegions::StdPrismExp::v_GetNfaces ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 787 of file StdPrismExp.cpp.

788  {
789  return 5;
790  }
int Nektar::StdRegions::StdPrismExp::v_GetNverts ( ) const
protectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 777 of file StdPrismExp.cpp.

778  {
779  return 6;
780  }
int Nektar::StdRegions::StdPrismExp::v_GetTotalEdgeIntNcoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 860 of file StdPrismExp.cpp.

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

861  {
862  int P = GetBasisNumModes(0)-2;
863  int Q = GetBasisNumModes(1)-2;
864  int R = GetBasisNumModes(2)-2;
865 
866  return 2*P+3*Q+3*R;
867  }
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:178
int Nektar::StdRegions::StdPrismExp::v_GetTotalFaceIntNcoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 909 of file StdPrismExp.cpp.

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

910  {
911  int Pi = GetBasisNumModes(0) - 2;
912  int Qi = GetBasisNumModes(1) - 2;
913  int Ri = GetBasisNumModes(2) - 2;
914 
915  return Pi * Qi +
916  Pi * (2*Ri - Pi - 1) +
917  2* Qi * Ri;
918  }
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:178
int Nektar::StdRegions::StdPrismExp::v_GetVertexMap ( int  localVertexId,
bool  useCoeffPacking = false 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalPrismExp.

Definition at line 1354 of file StdPrismExp.cpp.

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

1355  {
1359  "Mapping not defined for this type of basis");
1360 
1361  int l = 0;
1362 
1363  if(useCoeffPacking == true) // follow packing of coefficients i.e q,r,p
1364  {
1365  switch (vId)
1366  {
1367  case 0:
1368  l = GetMode(0,0,0);
1369  break;
1370  case 1:
1371  l = GetMode(0,0,1);
1372  break;
1373  case 2:
1374  l = GetMode(0,1,0);
1375  break;
1376  case 3:
1377  l = GetMode(0,1,1);
1378  break;
1379  case 4:
1380  l = GetMode(1,0,0);
1381  break;
1382  case 5:
1383  l = GetMode(1,1,0);
1384  break;
1385  default:
1386  ASSERTL0(false, "local vertex id must be between 0 and 5");
1387  }
1388  }
1389  else
1390  {
1391  switch (vId)
1392  {
1393  case 0:
1394  l = GetMode(0,0,0);
1395  break;
1396  case 1:
1397  l = GetMode(1,0,0);
1398  break;
1399  case 2:
1400  l = GetMode(1,1,0);
1401  break;
1402  case 3:
1403  l = GetMode(0,1,0);
1404  break;
1405  case 4:
1406  l = GetMode(0,0,1);
1407  break;
1408  case 5:
1409  l = GetMode(0,1,1);
1410  break;
1411  default:
1412  ASSERTL0(false, "local vertex id must be between 0 and 5");
1413  }
1414  }
1415 
1416  return l;
1417  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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
Principle Modified Functions .
Definition: BasisType.h:50
int GetMode(int I, int J, int K)
Compute the local mode number in the expansion for a particular tensorial combination.
void Nektar::StdRegions::StdPrismExp::v_IProductWRTBase ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

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

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

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

Implements Nektar::StdRegions::StdExpansion.

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

Definition at line 397 of file StdPrismExp.cpp.

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

Referenced by v_FwdTrans().

400  {
403  "Basis[1] is not a general tensor type");
404 
407  "Basis[2] is not a general tensor type");
408 
409  if(m_base[0]->Collocation() && m_base[1]->Collocation())
410  {
411  MultiplyByQuadratureMetric(inarray,outarray);
412  }
413  else
414  {
415  StdPrismExp::v_IProductWRTBase_SumFac(inarray,outarray);
416  }
417  }
Principle Modified Functions .
Definition: BasisType.h:51
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:947
Principle Orthogonal Functions .
Definition: BasisType.h:47
Principle Modified Functions .
Definition: BasisType.h:50
Principle Orthogonal Functions .
Definition: BasisType.h:48
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdPrismExp::v_IProductWRTBase_MatOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Implementation of the local matrix inner product operation.

Definition at line 422 of file StdPrismExp.cpp.

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

425  {
426  int nq = GetTotPoints();
427  StdMatrixKey iprodmatkey(eIProductWRTBase,DetShapeType(),*this);
428  DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
429 
430  Blas::Dgemv('N',m_ncoeffs,nq,1.0,iprodmat->GetPtr().get(),
431  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
432  }
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:706
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:141
void Nektar::StdRegions::StdPrismExp::v_IProductWRTBase_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
bool  multiplybyweights = true 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 434 of file StdPrismExp.cpp.

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

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

438  {
439  int nquad1 = m_base[1]->GetNumPoints();
440  int nquad2 = m_base[2]->GetNumPoints();
441  int order0 = m_base[0]->GetNumModes();
442  int order1 = m_base[1]->GetNumModes();
443 
444  Array<OneD, NekDouble> wsp(order0*nquad2*(nquad1+order1));
445 
446  if(multiplybyweights)
447  {
448  Array<OneD, NekDouble> tmp(inarray.num_elements());
449 
450  MultiplyByQuadratureMetric(inarray,tmp);
451  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
452  m_base[1]->GetBdata(),
453  m_base[2]->GetBdata(),
454  tmp,outarray,wsp,
455  true,true,true);
456  }
457  else
458  {
459  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
460  m_base[1]->GetBdata(),
461  m_base[2]->GetBdata(),
462  inarray,outarray,wsp,
463  true,true,true);
464  }
465  }
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:947
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)
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdPrismExp::v_IProductWRTBase_SumFacKernel ( const Array< OneD, const NekDouble > &  base0,
const Array< OneD, const NekDouble > &  base1,
const Array< OneD, const NekDouble > &  base2,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
Array< OneD, NekDouble > &  wsp,
bool  doCheckCollDir0,
bool  doCheckCollDir1,
bool  doCheckCollDir2 
)
protectedvirtual

Implements Nektar::StdRegions::StdExpansion3D.

Definition at line 467 of file StdPrismExp.cpp.

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

477  {
478  // Interior prism implementation based on Spen's book page
479  // 119. and 608.
480  const int nquad0 = m_base[0]->GetNumPoints();
481  const int nquad1 = m_base[1]->GetNumPoints();
482  const int nquad2 = m_base[2]->GetNumPoints();
483  const int order0 = m_base[0]->GetNumModes ();
484  const int order1 = m_base[1]->GetNumModes ();
485  const int order2 = m_base[2]->GetNumModes ();
486 
487  int i, mode;
488 
489  ASSERTL1(wsp.num_elements() >= nquad1*nquad2*order0 +
490  nquad2*order0*order1,
491  "Insufficient workspace size");
492 
493  Array<OneD, NekDouble> tmp0 = wsp;
494  Array<OneD, NekDouble> tmp1 = wsp + nquad1*nquad2*order0;
495 
496  // Inner product with respect to the '0' direction
497  Blas::Dgemm('T', 'N', nquad1*nquad2, order0, nquad0,
498  1.0, inarray.get(), nquad0,
499  base0.get(), nquad0,
500  0.0, tmp0.get(), nquad1*nquad2);
501 
502  // Inner product with respect to the '1' direction
503  Blas::Dgemm('T', 'N', nquad2*order0, order1, nquad1,
504  1.0, tmp0.get(), nquad1,
505  base1.get(), nquad1,
506  0.0, tmp1.get(), nquad2*order0);
507 
508  // Inner product with respect to the '2' direction
509  for (mode=i=0; i < order0; ++i)
510  {
511  Blas::Dgemm('T', 'N', order2-i, order1, nquad2,
512  1.0, base2.get() + mode*nquad2, nquad2,
513  tmp1.get() + i*nquad2, nquad2*order0,
514  0.0, outarray.get()+mode*order1, order2-i);
515  mode += order2-i;
516  }
517 
518  // Fix top singular vertices; performs phi_{0,q,1} +=
519  // phi_1(xi_1)*phi_q(xi_2)*phi_{01}*phi_r(xi_2).
521  {
522  for (i = 0; i < order1; ++i)
523  {
524  mode = GetMode(0,i,1);
525  outarray[mode] += Blas::Ddot(
526  nquad2, base2.get()+nquad2, 1,
527  tmp1.get()+i*order0*nquad2+nquad2, 1);
528  }
529  }
530  }
Principle Modified Functions .
Definition: BasisType.h:49
T Ddot(int n, const Array< OneD, const T > &w, const int incw, const Array< OneD, const T > &x, const int incx, const Array< OneD, const int > &y, const int incy)
Definition: VmathArray.hpp:436
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
int GetMode(int I, int J, int K)
Compute the local mode number in the expansion for a particular tensorial combination.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdPrismExp::v_IProductWRTDerivBase ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 536 of file StdPrismExp.cpp.

References v_IProductWRTDerivBase_SumFac().

540  {
541  v_IProductWRTDerivBase_SumFac(dir,inarray,outarray);
542  }
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void Nektar::StdRegions::StdPrismExp::v_IProductWRTDerivBase_MatOp ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Definition at line 544 of file StdPrismExp.cpp.

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

548  {
549  ASSERTL0(dir >= 0 && dir <= 2, "input dir is out of range");
550 
551  int nq = GetTotPoints();
552  MatrixType mtype;
553 
554  switch (dir)
555  {
556  case 0:
557  mtype = eIProductWRTDerivBase0;
558  break;
559  case 1:
560  mtype = eIProductWRTDerivBase1;
561  break;
562  case 2:
563  mtype = eIProductWRTDerivBase2;
564  break;
565  }
566 
567  StdMatrixKey iprodmatkey(mtype,DetShapeType(),*this);
568  DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
569 
570  Blas::Dgemv('N',m_ncoeffs,nq,1.0,iprodmat->GetPtr().get(),
571  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
572  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:470
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:706
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:141
void Nektar::StdRegions::StdPrismExp::v_IProductWRTDerivBase_SumFac ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 574 of file StdPrismExp.cpp.

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

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

578  {
579  ASSERTL0(dir >= 0 && dir <= 2, "input dir is out of range");
580 
581  int i;
582  int order0 = m_base[0]->GetNumModes ();
583  int order1 = m_base[1]->GetNumModes ();
584  int nquad0 = m_base[0]->GetNumPoints();
585  int nquad1 = m_base[1]->GetNumPoints();
586  int nquad2 = m_base[2]->GetNumPoints();
587 
588  const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
589  const Array<OneD, const NekDouble> &z2 = m_base[2]->GetZ();
590  Array<OneD, NekDouble> gfac0(nquad0);
591  Array<OneD, NekDouble> gfac2(nquad2);
592  Array<OneD, NekDouble> tmp0 (nquad0*nquad1*nquad2);
593  Array<OneD, NekDouble> wsp (order0*nquad2*(nquad1+order1));
594 
595  // set up geometric factor: (1+z0)/2
596  for (i = 0; i < nquad0; ++i)
597  {
598  gfac0[i] = 0.5*(1+z0[i]);
599  }
600 
601  // Set up geometric factor: 2/(1-z2)
602  for (i = 0; i < nquad2; ++i)
603  {
604  gfac2[i] = 2.0/(1-z2[i]);
605  }
606 
607  // Scale first derivative term by gfac2.
608  if (dir != 1)
609  {
610  for (i = 0; i < nquad2; ++i)
611  {
612  Vmath::Smul(nquad0*nquad1,gfac2[i],
613  &inarray[0]+i*nquad0*nquad1,1,
614  &tmp0 [0]+i*nquad0*nquad1,1);
615  }
616  MultiplyByQuadratureMetric(tmp0,tmp0);
617  }
618 
619  switch (dir)
620  {
621  case 0:
622  {
623  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
624  m_base[1]->GetBdata (),
625  m_base[2]->GetBdata (),
626  tmp0,outarray,wsp,
627  true,true,true);
628  break;
629  }
630 
631  case 1:
632  {
633  MultiplyByQuadratureMetric(inarray,tmp0);
634  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata (),
635  m_base[1]->GetDbdata(),
636  m_base[2]->GetBdata (),
637  tmp0,outarray,wsp,
638  true,true,true);
639  break;
640  }
641 
642  case 2:
643  {
644  Array<OneD, NekDouble> tmp1(m_ncoeffs);
645 
646  // Scale eta_1 derivative with gfac0.
647  for(i = 0; i < nquad1*nquad2; ++i)
648  {
649  Vmath::Vmul(nquad0,&gfac0[0],1,&tmp0[0]+i*nquad0,1,&tmp0[0]+i*nquad0,1);
650  }
651 
652  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
653  m_base[1]->GetBdata(),
654  m_base[2]->GetBdata(),
655  tmp0,tmp1,wsp,
656  true,true,true);
657 
658  MultiplyByQuadratureMetric(inarray, tmp0);
659  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
660  m_base[1]->GetBdata(),
661  m_base[2]->GetDbdata(),
662  tmp0,outarray,wsp,
663  true,true,true);
664 
665  Vmath::Vadd(m_ncoeffs,&tmp1[0],1,&outarray[0],1,&outarray[0],1);
666  break;
667  }
668  }
669  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:947
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:213
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:299
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:183
bool Nektar::StdRegions::StdPrismExp::v_IsBoundaryInteriorExpansion ( )
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1029 of file StdPrismExp.cpp.

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

1030  {
1031  return (m_base[0]->GetBasisType() == LibUtilities::eModified_A) &&
1032  (m_base[1]->GetBasisType() == LibUtilities::eModified_A) &&
1034  }
Principle Modified Functions .
Definition: BasisType.h:49
Principle Modified Functions .
Definition: BasisType.h:50
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdPrismExp::v_LocCoordToLocCollapsed ( const Array< OneD, const NekDouble > &  xi,
Array< OneD, NekDouble > &  eta 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 678 of file StdPrismExp.cpp.

References Nektar::NekConstants::kNekZeroTol.

681  {
682 
683  if( fabs(xi[2]-1.0) < NekConstants::kNekZeroTol)
684  {
685  // Very top point of the prism
686  eta[0] = -1.0;
687  eta[1] = xi[1];
688  eta[2] = 1.0;
689  }
690  else
691  {
692  // Third basis function collapsed to "pr" direction instead of
693  // "qr" direction
694  eta[2] = xi[2]; // eta_z = xi_z
695  eta[1] = xi[1]; //eta_y = xi_y
696  eta[0] = 2.0*(1.0 + xi[0])/(1.0 - xi[2]) - 1.0;
697  }
698  }
static const NekDouble kNekZeroTol
void Nektar::StdRegions::StdPrismExp::v_MultiplyByStdQuadratureMetric ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1942 of file StdPrismExp.cpp.

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

1945  {
1946  int i, j;
1947  int nquad0 = m_base[0]->GetNumPoints();
1948  int nquad1 = m_base[1]->GetNumPoints();
1949  int nquad2 = m_base[2]->GetNumPoints();
1950 
1951  const Array<OneD, const NekDouble>& w0 = m_base[0]->GetW();
1952  const Array<OneD, const NekDouble>& w1 = m_base[1]->GetW();
1953  const Array<OneD, const NekDouble>& w2 = m_base[2]->GetW();
1954 
1955  const Array<OneD, const NekDouble>& z2 = m_base[2]->GetZ();
1956 
1957  // Multiply by integration constants in x-direction
1958  for(i = 0; i < nquad1*nquad2; ++i)
1959  {
1960  Vmath::Vmul(nquad0, inarray.get()+i*nquad0, 1,
1961  w0.get(), 1, outarray.get()+i*nquad0,1);
1962  }
1963 
1964  // Multiply by integration constants in y-direction
1965  for(j = 0; j < nquad2; ++j)
1966  {
1967  for(i = 0; i < nquad1; ++i)
1968  {
1969  Blas::Dscal(nquad0,w1[i], &outarray[0]+i*nquad0 +
1970  j*nquad0*nquad1,1);
1971  }
1972  }
1973 
1974  // Multiply by integration constants in z-direction; need to
1975  // incorporate factor (1-eta_3)/2 into weights, but only if using
1976  // GLL quadrature points.
1977  switch(m_base[2]->GetPointsType())
1978  {
1979  // (1,0) Jacobi inner product.
1981  for(i = 0; i < nquad2; ++i)
1982  {
1983  Blas::Dscal(nquad0*nquad1, 0.5*w2[i],
1984  &outarray[0]+i*nquad0*nquad1, 1);
1985  }
1986  break;
1987 
1988  default:
1989  for(i = 0; i < nquad2; ++i)
1990  {
1991  Blas::Dscal(nquad0*nquad1,0.25*(1-z2[i])*w2[i],
1992  &outarray[0]+i*nquad0*nquad1,1);
1993  }
1994  break;
1995  }
1996 
1997  }
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:59
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:183
int Nektar::StdRegions::StdPrismExp::v_NumBndryCoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 801 of file StdPrismExp.cpp.

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

802  {
805  "BasisType is not a boundary interior form");
808  "BasisType is not a boundary interior form");
811  "BasisType is not a boundary interior form");
812 
813  int P = m_base[0]->GetNumModes();
814  int Q = m_base[1]->GetNumModes();
815  int R = m_base[2]->GetNumModes();
816 
819  }
Principle Modified Functions .
Definition: BasisType.h:49
Principle Modified Functions .
Definition: BasisType.h:50
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
int getNumberOfBndCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:297
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:228
Array< OneD, LibUtilities::BasisSharedPtr > m_base
int Nektar::StdRegions::StdPrismExp::v_NumDGBndryCoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 821 of file StdPrismExp.cpp.

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

822  {
825  "BasisType is not a boundary interior form");
828  "BasisType is not a boundary interior form");
831  "BasisType is not a boundary interior form");
832 
833  int P = m_base[0]->GetNumModes()-1;
834  int Q = m_base[1]->GetNumModes()-1;
835  int R = m_base[2]->GetNumModes()-1;
836 
837  return (P+1)*(Q+1) // 1 rect. face on base
838  + 2*(Q+1)*(R+1) // other 2 rect. faces
839  + 2*(R+1) + P*(1 + 2*R - P); // 2 tri. faces
840  }
Principle Modified Functions .
Definition: BasisType.h:49
Principle Modified Functions .
Definition: BasisType.h:50
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:228
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdPrismExp::v_PhysDeriv ( const Array< OneD, const NekDouble > &  u_physical,
Array< OneD, NekDouble > &  out_dxi1,
Array< OneD, NekDouble > &  out_dxi2,
Array< OneD, NekDouble > &  out_dxi3 
)
protectedvirtual

Calculate the derivative of the physical points.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PrismExp.

Definition at line 98 of file StdPrismExp.cpp.

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

Referenced by v_PhysDeriv(), and v_StdPhysDeriv().

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

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

See also
StdRegions::StdExpansion::PhysDeriv

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 162 of file StdPrismExp.cpp.

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

165  {
166  switch(dir)
167  {
168  case 0:
169  {
170  v_PhysDeriv(inarray, outarray, NullNekDouble1DArray,
172  break;
173  }
174 
175  case 1:
176  {
177  v_PhysDeriv(inarray, NullNekDouble1DArray, outarray,
179  break;
180  }
181 
182  case 2:
183  {
185  NullNekDouble1DArray, outarray);
186  break;
187  }
188 
189  default:
190  {
191  ASSERTL1(false,"input dir is out of range");
192  }
193  break;
194  }
195  }
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: StdPrismExp.cpp:98
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
void Nektar::StdRegions::StdPrismExp::v_ReduceOrderCoeffs ( int  numMin,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2066 of file StdPrismExp.cpp.

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

2070  {
2071  int nquad0 = m_base[0]->GetNumPoints();
2072  int nquad1 = m_base[1]->GetNumPoints();
2073  int nquad2 = m_base[2]->GetNumPoints();
2074  int nqtot = nquad0*nquad1*nquad2;
2075  int nmodes0 = m_base[0]->GetNumModes();
2076  int nmodes1 = m_base[1]->GetNumModes();
2077  int nmodes2 = m_base[2]->GetNumModes();
2078  int numMax = nmodes0;
2079 
2080  Array<OneD, NekDouble> coeff (m_ncoeffs);
2081  Array<OneD, NekDouble> coeff_tmp1(m_ncoeffs, 0.0);
2082  Array<OneD, NekDouble> phys_tmp (nqtot, 0.0);
2083  Array<OneD, NekDouble> tmp, tmp2, tmp3, tmp4;
2084 
2085 
2086  const LibUtilities::PointsKey Pkey0 = m_base[0]->GetPointsKey();
2087  const LibUtilities::PointsKey Pkey1 = m_base[1]->GetPointsKey();
2088  const LibUtilities::PointsKey Pkey2 = m_base[2]->GetPointsKey();
2089 
2090  LibUtilities::BasisKey bortho0(
2091  LibUtilities::eOrtho_A, nmodes0, Pkey0);
2092  LibUtilities::BasisKey bortho1(
2093  LibUtilities::eOrtho_A, nmodes1, Pkey1);
2094  LibUtilities::BasisKey bortho2(
2095  LibUtilities::eOrtho_B, nmodes2, Pkey2);
2096 
2097  int cnt = 0;
2098  int u = 0;
2099  int i = 0;
2100  StdRegions::StdPrismExpSharedPtr OrthoPrismExp;
2101 
2103  ::AllocateSharedPtr(bortho0, bortho1, bortho2);
2104 
2105  BwdTrans(inarray,phys_tmp);
2106  OrthoPrismExp->FwdTrans(phys_tmp, coeff);
2107 
2108  // filtering
2109  for (u = 0; u < numMin; ++u)
2110  {
2111  for (i = 0; i < numMin; ++i)
2112  {
2113  Vmath::Vcopy(numMin - u, tmp = coeff + cnt, 1,
2114  tmp2 = coeff_tmp1 + cnt, 1);
2115  cnt += numMax - u;
2116  }
2117 
2118  for (i = numMin; i < numMax; ++i)
2119  {
2120  cnt += numMax - u;
2121  }
2122  }
2123 
2124  OrthoPrismExp->BwdTrans(coeff_tmp1, phys_tmp);
2125  StdPrismExp::FwdTrans(phys_tmp, outarray);
2126  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< StdPrismExp > StdPrismExpSharedPtr
Definition: StdPrismExp.h:258
Principle Orthogonal Functions .
Definition: BasisType.h:47
Principle Orthogonal Functions .
Definition: BasisType.h:46
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Backward transformation from coefficient space to physical space...
Definition: StdExpansion.h:531
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Forward transformation from physical space to coefficient space...
void Nektar::StdRegions::StdPrismExp::v_StdPhysDeriv ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d0,
Array< OneD, NekDouble > &  out_d1,
Array< OneD, NekDouble > &  out_d2 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 197 of file StdPrismExp.cpp.

References v_PhysDeriv().

201  {
202  StdPrismExp::v_PhysDeriv(inarray, out_d0, out_d1, out_d2);
203  }
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
Calculate the derivative of the physical points.
Definition: StdPrismExp.cpp:98
void Nektar::StdRegions::StdPrismExp::v_StdPhysDeriv ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 205 of file StdPrismExp.cpp.

References v_PhysDeriv().

208  {
209  StdPrismExp::v_PhysDeriv(dir, inarray, outarray);
210  }
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
Calculate the derivative of the physical points.
Definition: StdPrismExp.cpp:98
void Nektar::StdRegions::StdPrismExp::v_SVVLaplacianFilter ( Array< OneD, NekDouble > &  array,
const StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PrismExp.

Definition at line 1999 of file StdPrismExp.cpp.

References Nektar::StdRegions::StdExpansion::BwdTrans(), Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::eFactorSVVDiffCoeff, Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::eOrtho_B, Nektar::StdRegions::StdExpansion::FwdTrans(), Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), Nektar::StdRegions::StdExpansion::GetPointsType(), and Nektar::StdRegions::StdExpansion::m_base.

2001  {
2002  // Generate an orthonogal expansion
2003  int qa = m_base[0]->GetNumPoints();
2004  int qb = m_base[1]->GetNumPoints();
2005  int qc = m_base[2]->GetNumPoints();
2006  int nmodes_a = m_base[0]->GetNumModes();
2007  int nmodes_b = m_base[1]->GetNumModes();
2008  int nmodes_c = m_base[2]->GetNumModes();
2009  // Declare orthogonal basis.
2010  LibUtilities::PointsKey pa(qa,m_base[0]->GetPointsType());
2011  LibUtilities::PointsKey pb(qb,m_base[1]->GetPointsType());
2012  LibUtilities::PointsKey pc(qc,m_base[2]->GetPointsType());
2013 
2014  LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A,nmodes_a,pa);
2015  LibUtilities::BasisKey Bb(LibUtilities::eOrtho_A,nmodes_b,pb);
2016  LibUtilities::BasisKey Bc(LibUtilities::eOrtho_B,nmodes_c,pc);
2017  StdPrismExp OrthoExp(Ba,Bb,Bc);
2018 
2019  Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
2020  int i,j,k,cnt = 0;
2021 
2022  //SVV filter paramaters (how much added diffusion relative to physical one
2023  // and fraction of modes from which you start applying this added diffusion)
2024  //
2025  NekDouble SvvDiffCoeff = mkey.GetConstFactor(StdRegions::eFactorSVVDiffCoeff);
2026  NekDouble SVVCutOff = mkey.GetConstFactor(StdRegions::eFactorSVVCutoffRatio);
2027 
2028  //Defining the cut of mode
2029  int cutoff_a = (int) (SVVCutOff*nmodes_a);
2030  int cutoff_b = (int) (SVVCutOff*nmodes_b);
2031  int cutoff_c = (int) (SVVCutOff*nmodes_c);
2032  //To avoid the fac[j] from blowing up
2033  NekDouble epsilon = 1;
2034 
2035  // project onto modal space.
2036  OrthoExp.FwdTrans(array,orthocoeffs);
2037  int nmodes = min(min(nmodes_a,nmodes_b),nmodes_c);
2038  NekDouble cutoff = min(min(cutoff_a,cutoff_b),cutoff_c);
2039 
2040  //------"New" Version August 22nd '13--------------------
2041  for(i = 0; i < nmodes_a; ++i)//P
2042  {
2043  for(j = 0; j < nmodes_b; ++j) //Q
2044  {
2045  for(k = 0; k < nmodes_c-i; ++k) //R
2046  {
2047  if(j >= cutoff || i + k >= cutoff)
2048  {
2049  orthocoeffs[cnt] *= (SvvDiffCoeff*exp(-(i+k-nmodes)*(i+k-nmodes)/((NekDouble)((i+k-cutoff+epsilon)*(i+k-cutoff+epsilon))))*exp(-(j-nmodes)*(j-nmodes)/((NekDouble)((j-cutoff+epsilon)*(j-cutoff+epsilon)))));
2050  }
2051  else
2052  {
2053  orthocoeffs[cnt] *= 0.0;
2054  }
2055  cnt++;
2056  }
2057  }
2058  }
2059 
2060  // backward transform to physical space
2061  OrthoExp.BwdTrans(orthocoeffs,array);
2062  }
Principle Orthogonal Functions .
Definition: BasisType.h:47
Principle Orthogonal Functions .
Definition: BasisType.h:46
double NekDouble
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