Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator 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.
void BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
void IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion
 StdExpansion ()
 Default Constructor.
 StdExpansion (const int numcoeffs, const int numbases, const LibUtilities::BasisKey &Ba=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bb=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bc=LibUtilities::NullBasisKey)
 Constructor.
 StdExpansion (const StdExpansion &T)
 Copy Constructor.
virtual ~StdExpansion ()
 Destructor.
int GetNumBases () const
 This function returns the number of 1D bases used in the expansion.
const Array< OneD, const
LibUtilities::BasisSharedPtr > & 
GetBase () const
 This function gets the shared point to basis.
const
LibUtilities::BasisSharedPtr
GetBasis (int dir) const
 This function gets the shared point to basis in the dir direction.
int GetNcoeffs (void) const
 This function returns the total number of coefficients used in the expansion.
int GetTotPoints () const
 This function returns the total number of quadrature points used in the element.
LibUtilities::BasisType GetBasisType (const int dir) const
 This function returns the type of basis used in the dir direction.
int GetBasisNumModes (const int dir) const
 This function returns the number of expansion modes in the dir direction.
int EvalBasisNumModesMax (void) const
 This function returns the maximum number of expansion modes over all local directions.
LibUtilities::PointsType GetPointsType (const int dir) const
 This function returns the type of quadrature points used in the dir direction.
int GetNumPoints (const int dir) const
 This function returns the number of quadrature points in the dir direction.
const Array< OneD, const
NekDouble > & 
GetPoints (const int dir) const
 This function returns a pointer to the array containing the quadrature points in dir direction.
int GetNverts () const
 This function returns the number of vertices of the expansion domain.
int GetNedges () const
 This function returns the number of edges of the expansion domain.
int GetEdgeNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th edge.
int GetTotalEdgeIntNcoeffs () const
int GetEdgeNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th edge.
int DetCartesianDirOfEdge (const int edge)
const LibUtilities::BasisKey DetEdgeBasisKey (const int i) const
const LibUtilities::BasisKey DetFaceBasisKey (const int i, const int k) const
int GetFaceNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th face.
int GetFaceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th face.
int GetFaceIntNcoeffs (const int i) const
int GetTotalFaceIntNcoeffs () const
LibUtilities::PointsKey GetFacePointsKey (const int i, const int j) const
int NumBndryCoeffs (void) const
int NumDGBndryCoeffs (void) const
LibUtilities::BasisType GetEdgeBasisType (const int i) const
 This function returns the type of expansion basis on the i-th edge.
int GetNfaces () const
 This function returns the number of faces of the expansion domain.
int GetNtrace () const
 Returns the number of trace elements connected to this element.
LibUtilities::ShapeType DetShapeType () const
 This function returns the shape of the expansion domain.
int GetShapeDimension () const
bool IsBoundaryInteriorExpansion ()
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space.
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space.
void FwdTrans_BndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
NekDouble Integral (const Array< OneD, const NekDouble > &inarray)
 This function integrates the specified function over the domain.
void FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 This function fills the array outarray with the mode-th mode of the expansion.
void IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 this function calculates the inner product of a given function f with the different modes of the expansion
void IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check)
void IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
int GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id.
void SetElmtId (const int id)
 Set the element id of this expansion when used in a list by returning value of m_elmt_id.
void GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2=NullNekDouble1DArray, Array< OneD, NekDouble > &coords_3=NullNekDouble1DArray)
 this function returns the physical coordinates of the quadrature points of the expansion
void GetCoord (const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
 given the coordinates of a point of the element in the local collapsed coordinate system, this function calculates the physical coordinates of the point
DNekMatSharedPtr GetStdMatrix (const StdMatrixKey &mkey)
DNekBlkMatSharedPtr GetStdStaticCondMatrix (const StdMatrixKey &mkey)
IndexMapValuesSharedPtr GetIndexMap (const IndexMapKey &ikey)
const Array< OneD, const
NekDouble > & 
GetPhysNormals (void)
void SetPhysNormals (Array< OneD, const NekDouble > &normal)
virtual void SetUpPhysNormals (const int edge)
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
void DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
StdRegions::Orientation GetForient (int face)
StdRegions::Orientation GetEorient (int edge)
StdRegions::Orientation GetPorient (int point)
StdRegions::Orientation GetCartesianEorient (int edge)
void SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
void SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
int CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
void ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs)
NekDouble StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
int GetCoordim ()
void GetBoundaryMap (Array< OneD, unsigned int > &outarray)
void GetInteriorMap (Array< OneD, unsigned int > &outarray)
int GetVertexMap (const int localVertexId, bool useCoeffPacking=false)
void GetEdgeInteriorMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
void GetFaceInteriorMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
void GetEdgeToElementMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
void GetFaceToElementMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int nummodesA=-1, int nummodesB=-1)
void GetEdgePhysVals (const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Extract the physical values along edge edge from inarray into outarray following the local edge orientation and point distribution defined by defined in EdgeExp.
void GetEdgePhysVals (const int edge, const boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void GetTracePhysVals (const int edge, const boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void GetVertexPhysVals (const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray)
void GetEdgeInterpVals (const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void GetEdgeQFactors (const int edge, Array< OneD, NekDouble > &outarray)
 Extract the metric factors to compute the contravariant fluxes along edge edge and stores them into outarray following the local edge orientation (i.e. anticlockwise convention).
void GetFacePhysVals (const int face, const boost::shared_ptr< StdExpansion > &FaceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=eNoOrientation)
void MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
DNekMatSharedPtr CreateGeneralMatrix (const StdMatrixKey &mkey)
 this function generates the mass matrix $\mathbf{M}[i][j] = \int \phi_i(\mathbf{x}) \phi_j(\mathbf{x}) d\mathbf{x}$
void GeneralMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
void LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void WeakDirectionalDerivMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void MassLevelCurvatureMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LinearAdvectionDiffusionReactionMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
void HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
DNekMatSharedPtr GenMatrix (const StdMatrixKey &mkey)
void PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void PhysDeriv_s (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds)
void PhysDeriv_n (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn)
void PhysDirectionalDeriv (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &outarray)
void StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
void StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void AddRobinMassMatrix (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
void AddRobinEdgeContribution (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, Array< OneD, NekDouble > &coeffs)
NekDouble PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain.
NekDouble PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain.
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta.
const boost::shared_ptr
< SpatialDomains::GeomFactors > & 
GetMetricInfo (void) const
virtual int v_GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id.
virtual const Array< OneD,
const NekDouble > & 
v_GetPhysNormals (void)
virtual void v_SetPhysNormals (Array< OneD, const NekDouble > &normal)
virtual void v_SetUpPhysNormals (const int edge)
virtual void v_ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmode_offset, NekDouble *coeffs)
 Unpack data from input file assuming it comes from the same expansion type.
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
virtual void v_DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
virtual StdRegions::Orientation v_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.
NekDouble L2 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete $ L_2$ error, $ | \epsilon |_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 dx \right]^{1/2} d\xi_1 $ where $ u_{exact}$ is given by the array sol.
NekDouble H1 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete $ H^1$ error, $ | \epsilon |^1_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 + \nabla(u - u_{exact})\cdot\nabla(u - u_{exact})\cdot dx \right]^{1/2} d\xi_1 $ where $ u_{exact}$ is given by the array sol.
const NormalVectorGetEdgeNormal (const int edge) const
void ComputeEdgeNormal (const int edge)
void NegateEdgeNormal (const int edge)
bool EdgeNormalNegated (const int edge)
void ComputeFaceNormal (const int face)
void NegateFaceNormal (const int face)
void ComputeVertexNormal (const int vertex)
const NormalVectorGetFaceNormal (const int face) const
const NormalVectorGetVertexNormal (const int vertex) const
const NormalVectorGetSurfaceNormal (const int id) const
const LibUtilities::PointsKeyVector GetPointsKeys () const
Array< OneD, unsigned int > GetEdgeInverseBoundaryMap (int eid)
Array< OneD, unsigned int > GetFaceInverseBoundaryMap (int fid, StdRegions::Orientation faceOrient=eNoOrientation)
DNekMatSharedPtr BuildInverseTransformationMatrix (const DNekScalMatSharedPtr &m_transformationmatrix)
void PhysInterpToSimplexEquiSpaced (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 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.
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.
template<class T >
boost::shared_ptr< T > as ()

Protected Member Functions

virtual void v_PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
 Calculate the derivative of the physical points.
virtual void v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Calculate the derivative of the physical points in a given direction.
virtual void v_StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
virtual void v_StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
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.
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:
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)
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.
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 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.
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 nummodesA=-1, int nummodesB=-1)
virtual int v_GetVertexMap (int localVertexId, bool useCoeffPacking=false)
virtual void v_GetEdgeInteriorMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
virtual void v_GetFaceInteriorMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
virtual void v_GetInteriorMap (Array< OneD, unsigned int > &outarray)
virtual void v_GetBoundaryMap (Array< OneD, unsigned int > &outarray)
virtual DNekMatSharedPtr v_GenMatrix (const StdMatrixKey &mkey)
virtual DNekMatSharedPtr v_CreateStdMatrix (const StdMatrixKey &mkey)
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.
virtual NekDouble v_PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
virtual void v_LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
virtual void v_HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
virtual NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray)
 Integrates the specified function over the domain.
virtual void v_NegateFaceNormal (const int face)
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion
DNekMatSharedPtr CreateStdMatrix (const StdMatrixKey &mkey)
DNekBlkMatSharedPtr CreateStdStaticCondMatrix (const StdMatrixKey &mkey)
 Create the static condensation of a matrix when using a boundary interior decomposition.
IndexMapValuesSharedPtr CreateIndexMap (const IndexMapKey &ikey)
 Create an IndexMap which contains mapping information linking any specific element shape with either its boundaries, edges, faces, verteces, etc.
void BwdTrans_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void GeneralMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void MassMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
void LaplacianMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LaplacianMatrixOp_MatFree (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void WeakDerivMatrixOp_MatFree (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void WeakDirectionalDerivMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void MassLevelCurvatureMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LinearAdvectionDiffusionReactionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
void HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void HelmholtzMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
virtual NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)

Private Member Functions

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

Additional Inherited Members

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

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

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

Definition at line 48 of file StdPrismExp.cpp.

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

Ba.GetNumModes(),
Bb.GetNumModes(),
Bc.GetNumModes()),
3,Ba,Bb,Bc),
Ba.GetNumModes(),
Bb.GetNumModes(),
Bc.GetNumModes()),
Ba,Bb,Bc)
{
ASSERTL0(Ba.GetNumModes() <= Bc.GetNumModes(),
"order in 'a' direction is higher than order in 'c' direction");
}
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 66 of file StdPrismExp.cpp.

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

Definition at line 74 of file StdPrismExp.cpp.

{
}

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

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

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

{
int Q = m_base[1]->GetNumModes() - 1;
int R = m_base[2]->GetNumModes() - 1;
return r + // Skip along stacks (r-direction)
q*(R+1-p) + // Skip along columns (q-direction)
(Q+1)*(p*R + 1-(p-2)*(p-1)/2); // Skip along rows (p-direction)
}
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 240 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().

{
"Basis[1] is not a general tensor type");
"Basis[2] is not a general tensor type");
if(m_base[0]->Collocation() &&
m_base[1]->Collocation() &&
m_base[2]->Collocation())
{
inarray, 1, outarray, 1);
}
else
{
StdPrismExp::v_BwdTrans_SumFac(inarray,outarray);
}
}
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 266 of file StdPrismExp.cpp.

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

Referenced by v_BwdTrans().

{
int nquad1 = m_base[1]->GetNumPoints();
int nquad2 = m_base[2]->GetNumPoints();
int order0 = m_base[0]->GetNumModes();
int order1 = m_base[1]->GetNumModes();
Array<OneD, NekDouble> wsp(nquad2*order1*order0 +
nquad1*nquad2*order0);
BwdTrans_SumFacKernel(m_base[0]->GetBdata(),
m_base[1]->GetBdata(),
m_base[2]->GetBdata(),
inarray,outarray,wsp,true,true,true);
}
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 284 of file StdPrismExp.cpp.

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

{
int i, mode;
int nquad0 = m_base[0]->GetNumPoints();
int nquad1 = m_base[1]->GetNumPoints();
int nquad2 = m_base[2]->GetNumPoints();
int nummodes0 = m_base[0]->GetNumModes();
int nummodes1 = m_base[1]->GetNumModes();
int nummodes2 = m_base[2]->GetNumModes();
Array<OneD, NekDouble> tmp0 = wsp;
Array<OneD, NekDouble> tmp1 = tmp0 + nquad2*nummodes1*nummodes0;
for (i = mode = 0; i < nummodes0; ++i)
{
Blas::Dgemm('N', 'N', nquad2, nummodes1, nummodes2-i,
1.0, base2.get() + mode*nquad2, nquad2,
inarray.get() + mode*nummodes1, nummodes2-i,
0.0, tmp0.get() + i*nquad2*nummodes1, nquad2);
mode += nummodes2-i;
}
{
for(i = 0; i < nummodes1; i++)
{
Blas::Daxpy(nquad2,inarray[1+i*nummodes2],base2.get()+nquad2,1,
tmp0.get()+nquad2*(nummodes1+i),1);
}
}
for (i = 0; i < nummodes0; i++)
{
Blas::Dgemm('N', 'T', nquad1, nquad2, nummodes1,
1.0, base1.get(), nquad1,
tmp0.get() + i*nquad2*nummodes1, nquad2,
0.0, tmp1.get() + i*nquad2*nquad1, nquad1);
}
Blas::Dgemm('N', 'T', nquad0, nquad2*nquad1, nummodes0,
1.0, base0.get(), nquad0,
tmp1.get(), nquad2*nquad1,
0.0, outarray.get(), nquad0);
}
int Nektar::StdRegions::StdPrismExp::v_CalcNumberOfCoefficients ( const std::vector< unsigned int > &  nummodes,
int &  modes_offset 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 945 of file StdPrismExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 1776 of file StdPrismExp.cpp.

References v_GenMatrix().

{
return v_GenMatrix(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 906 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().

{
ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
ASSERTL2(k >= 0 && k <= 1, "basis key id is out of range");
switch(i)
{
case 0:
{
m_base[k]->GetNumModes());
}
case 2:
case 4:
{
m_base[k+1]->GetNumModes());
}
case 1:
case 3:
{
m_base[2*k]->GetNumModes());
}
break;
}
// Should never get here.
}
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 741 of file StdPrismExp.cpp.

References Nektar::LibUtilities::ePrism.

{
}
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 710 of file StdPrismExp.cpp.

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

{
Array<OneD, NekDouble> tmp(m_ncoeffs,0.0);
tmp[mode] = 1.0;
StdPrismExp::v_BwdTrans(tmp, 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 348 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().

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 1699 of file StdPrismExp.cpp.

References Nektar::StdRegions::StdExpansion::CreateGeneralMatrix(), Nektar::StdRegions::ePhysInterpToEquiSpaced, Nektar::StdRegions::StdMatrixKey::GetMatrixType(), Nektar::LibUtilities::StdSegData::getNumberOfCoefficients(), Nektar::StdRegions::StdExpansion::LocCoordToLocCollapsed(), Nektar::StdRegions::StdExpansion::m_base, Vmath::Smul(), and Vmath::Vcopy().

Referenced by v_CreateStdMatrix().

{
MatrixType mtype = mkey.GetMatrixType();
switch(mtype)
{
{
int nq0 = m_base[0]->GetNumPoints();
int nq1 = m_base[1]->GetNumPoints();
int nq2 = m_base[2]->GetNumPoints();
int nq = max(nq0,max(nq1,nq2));
Array<OneD, Array<OneD, NekDouble> > coords (neq);
Array<OneD, NekDouble> coll (3);
Array<OneD, DNekMatSharedPtr> I (3);
Array<OneD, NekDouble> tmp (nq0);
Mat = MemoryManager<DNekMat>::
AllocateSharedPtr(neq,nq0*nq1*nq2);
int cnt = 0;
for(int i = 0; i < nq; ++i)
{
for(int j = 0; j < nq; ++j)
{
for(int k = 0; k < nq-i; ++k,++cnt)
{
coords[cnt] = Array<OneD, NekDouble>(3);
coords[cnt][0] = -1.0 + 2*k/(NekDouble)(nq-1);
coords[cnt][1] = -1.0 + 2*j/(NekDouble)(nq-1);
coords[cnt][2] = -1.0 + 2*i/(NekDouble)(nq-1);
}
}
}
for(int i = 0; i < neq; ++i)
{
LocCoordToLocCollapsed(coords[i],coll);
I[0] = m_base[0]->GetI(coll );
I[1] = m_base[1]->GetI(coll+1);
I[2] = m_base[2]->GetI(coll+2);
// interpolate first coordinate direction
NekDouble fac;
for( int k = 0; k < nq2; ++k)
{
for (int j = 0; j < nq1; ++j)
{
fac = (I[1]->GetPtr())[j]*(I[2]->GetPtr())[k];
Vmath::Smul(nq0,fac,I[0]->GetPtr(),1,tmp,1);
Vmath::Vcopy(nq0, &tmp[0], 1,
Mat->GetRawPtr() +
k * nq0 * nq1 * neq +
j * nq0 * neq + i,
neq);
}
}
}
}
break;
default:
{
}
break;
}
return Mat;
}
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 1640 of file StdPrismExp.cpp.

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

{
"BasisType is not a boundary interior form");
"BasisType is not a boundary interior form");
"BasisType is not a boundary interior form");
int P = m_base[0]->GetNumModes() - 1, p;
int Q = m_base[1]->GetNumModes() - 1, q;
int R = m_base[2]->GetNumModes() - 1, r;
int idx = 0;
// Loop over all boundary modes (in ascending order).
for (p = 0; p <= P; ++p)
{
// First two q-r planes are entirely boundary modes.
if (p <= 1)
{
for (q = 0; q <= Q; ++q)
{
for (r = 0; r <= R-p; ++r)
{
maparray[idx++] = GetMode(p,q,r);
}
}
}
else
{
// Remaining q-r planes contain boundary modes on the two
// left-hand sides and bottom edge.
for (q = 0; q <= Q; ++q)
{
if (q <= 1)
{
for (r = 0; r <= R-p; ++r)
{
maparray[idx++] = GetMode(p,q,r);
}
}
else
{
maparray[idx++] = GetMode(p,q,0);
}
}
}
}
}
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 686 of file StdPrismExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 957 of file StdPrismExp.cpp.

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

{
ASSERTL2(i >= 0 && i <= 8, "edge id is out of range");
if (i == 0 || i == 2)
{
return GetBasisType(0);
}
else if (i == 1 || i == 3 || i == 8)
{
return GetBasisType(1);
}
else
{
return GetBasisType(2);
}
}
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 1306 of file StdPrismExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 787 of file StdPrismExp.cpp.

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

Referenced by v_GetEdgeInteriorMap().

{
ASSERTL2(i >= 0 && i <= 8, "edge id is out of range");
if (i == 0 || i == 2)
{
return GetBasisNumModes(0);
}
else if (i == 1 || i == 3 || i == 8)
{
return GetBasisNumModes(1);
}
else
{
return GetBasisNumModes(2);
}
}
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 1420 of file StdPrismExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 832 of file StdPrismExp.cpp.

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

Referenced by v_GetFaceInteriorMap().

{
ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
int Pi = GetBasisNumModes(0) - 2;
int Qi = GetBasisNumModes(1) - 2;
int Ri = GetBasisNumModes(2) - 2;
if (i == 0)
{
return Pi * Qi;
}
else if (i == 1 || i == 3)
{
return Pi * (2*Ri - Pi - 1) / 2;
}
else
{
return Qi * Ri;
}
}
int Nektar::StdRegions::StdPrismExp::v_GetFaceNcoeffs ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 814 of file StdPrismExp.cpp.

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

{
ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
if (i == 0)
{
}
else if (i == 1 || i == 3)
{
int P = GetBasisNumModes(0)-1, Q = GetBasisNumModes(2)-1;
return Q+1 + (P*(1 + 2*Q - P))/2;
}
else
{
}
}
int Nektar::StdRegions::StdPrismExp::v_GetFaceNumPoints ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 865 of file StdPrismExp.cpp.

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

{
ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
if (i == 0)
{
return m_base[0]->GetNumPoints()*
m_base[1]->GetNumPoints();
}
else if (i == 1 || i == 3)
{
return m_base[0]->GetNumPoints()*
m_base[2]->GetNumPoints();
}
else
{
return m_base[1]->GetNumPoints()*
m_base[2]->GetNumPoints();
}
}
LibUtilities::PointsKey Nektar::StdRegions::StdPrismExp::v_GetFacePointsKey ( const int  i,
const int  j 
) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 886 of file StdPrismExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 986 of file StdPrismExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalPrismExp.

Definition at line 1602 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, and Nektar::StdRegions::StdExpansion::NumBndryCoeffs().

{
"BasisType is not a boundary interior form");
"BasisType is not a boundary interior form");
"BasisType is not a boundary interior form");
int P = m_base[0]->GetNumModes() - 1, p;
int Q = m_base[1]->GetNumModes() - 1, q;
int R = m_base[2]->GetNumModes() - 1, r;
int nIntCoeffs = m_ncoeffs - NumBndryCoeffs();
if(outarray.num_elements()!=nIntCoeffs)
{
outarray = Array<OneD, unsigned int>(nIntCoeffs);
}
int idx = 0;
// Loop over all interior modes.
for (p = 2; p <= P; ++p)
{
for (q = 2; q <= Q; ++q)
{
for (r = 1; r <= R-p; ++r)
{
outarray[idx++] = GetMode(p,q,r);
}
}
}
}
int Nektar::StdRegions::StdPrismExp::v_GetNedges ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 727 of file StdPrismExp.cpp.

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 732 of file StdPrismExp.cpp.

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 722 of file StdPrismExp.cpp.

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 805 of file StdPrismExp.cpp.

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

{
int P = GetBasisNumModes(0)-2;
int Q = GetBasisNumModes(1)-2;
int R = GetBasisNumModes(2)-2;
return 2*P+3*Q+3*R;
}
int Nektar::StdRegions::StdPrismExp::v_GetTotalFaceIntNcoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 854 of file StdPrismExp.cpp.

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

{
int Pi = GetBasisNumModes(0) - 2;
int Qi = GetBasisNumModes(1) - 2;
int Ri = GetBasisNumModes(2) - 2;
return Pi * Qi +
Pi * (2*Ri - Pi - 1) +
2* Qi * Ri;
}
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 1241 of file StdPrismExp.cpp.

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

{
"Mapping not defined for this type of basis");
int l = 0;
if(useCoeffPacking == true) // follow packing of coefficients i.e q,r,p
{
switch (vId)
{
case 0:
l = GetMode(0,0,0);
break;
case 1:
l = GetMode(0,0,1);
break;
case 2:
l = GetMode(0,1,0);
break;
case 3:
l = GetMode(0,1,1);
break;
case 4:
l = GetMode(1,0,0);
break;
case 5:
l = GetMode(1,1,0);
break;
default:
ASSERTL0(false, "local vertex id must be between 0 and 5");
}
}
else
{
switch (vId)
{
case 0:
l = GetMode(0,0,0);
break;
case 1:
l = GetMode(1,0,0);
break;
case 2:
l = GetMode(1,1,0);
break;
case 3:
l = GetMode(0,1,0);
break;
case 4:
l = GetMode(0,0,1);
break;
case 5:
l = GetMode(0,1,1);
break;
default:
ASSERTL0(false, "local vertex id must be between 0 and 5");
}
}
return l;
}
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 395 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().

{
"Basis[1] is not a general tensor type");
"Basis[2] is not a general tensor type");
if(m_base[0]->Collocation() && m_base[1]->Collocation())
{
MultiplyByQuadratureMetric(inarray,outarray);
}
else
{
}
}
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 420 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.

{
int nq = GetTotPoints();
StdMatrixKey iprodmatkey(eIProductWRTBase,DetShapeType(),*this);
DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
Blas::Dgemv('N',m_ncoeffs,nq,1.0,iprodmat->GetPtr().get(),
m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
}
void Nektar::StdRegions::StdPrismExp::v_IProductWRTBase_SumFac ( 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 432 of file StdPrismExp.cpp.

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

Referenced by v_IProductWRTBase().

{
int nquad1 = m_base[1]->GetNumPoints();
int nquad2 = m_base[2]->GetNumPoints();
int order0 = m_base[0]->GetNumModes();
int order1 = m_base[1]->GetNumModes();
Array<OneD, NekDouble> tmp(inarray.num_elements());
Array<OneD, NekDouble> wsp(order0*nquad2*(nquad1+order1));
m_base[1]->GetBdata(),
m_base[2]->GetBdata(),
tmp,outarray,wsp,
true,true,true);
}
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 453 of file StdPrismExp.cpp.

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

{
// Interior prism implementation based on Spen's book page
// 119. and 608.
const int nquad0 = m_base[0]->GetNumPoints();
const int nquad1 = m_base[1]->GetNumPoints();
const int nquad2 = m_base[2]->GetNumPoints();
const int order0 = m_base[0]->GetNumModes ();
const int order1 = m_base[1]->GetNumModes ();
const int order2 = m_base[2]->GetNumModes ();
int i, mode;
ASSERTL1(wsp.num_elements() >= nquad1*nquad2*order0 +
nquad2*order0*order1,
"Insufficient workspace size");
Array<OneD, NekDouble> tmp0 = wsp;
Array<OneD, NekDouble> tmp1 = wsp + nquad1*nquad2*order0;
// Inner product with respect to the '0' direction
Blas::Dgemm('T', 'N', nquad1*nquad2, order0, nquad0,
1.0, inarray.get(), nquad0,
base0.get(), nquad0,
0.0, tmp0.get(), nquad1*nquad2);
// Inner product with respect to the '1' direction
Blas::Dgemm('T', 'N', nquad2*order0, order1, nquad1,
1.0, tmp0.get(), nquad1,
base1.get(), nquad1,
0.0, tmp1.get(), nquad2*order0);
// Inner product with respect to the '2' direction
for (mode=i=0; i < order0; ++i)
{
Blas::Dgemm('T', 'N', order2-i, order1, nquad2,
1.0, base2.get() + mode*nquad2, nquad2,
tmp1.get() + i*nquad2, nquad2*order0,
0.0, outarray.get()+mode*order1, order2-i);
mode += order2-i;
}
// Fix top singular vertices; performs phi_{0,q,1} +=
// phi_1(xi_1)*phi_q(xi_2)*phi_{01}*phi_r(xi_2).
{
for (i = 0; i < order1; ++i)
{
mode = GetMode(0,i,1);
outarray[mode] += Blas::Ddot(
nquad2, base2.get()+nquad2, 1,
tmp1.get()+i*order0*nquad2+nquad2, 1);
}
}
}
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 522 of file StdPrismExp.cpp.

References v_IProductWRTDerivBase_SumFac().

{
v_IProductWRTDerivBase_SumFac(dir,inarray,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 530 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.

{
ASSERTL0(dir >= 0 && dir <= 2, "input dir is out of range");
int nq = GetTotPoints();
MatrixType mtype;
switch (dir)
{
case 0:
break;
case 1:
break;
case 2:
break;
}
StdMatrixKey iprodmatkey(mtype,DetShapeType(),*this);
DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
Blas::Dgemv('N',m_ncoeffs,nq,1.0,iprodmat->GetPtr().get(),
m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
}
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 560 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().

{
ASSERTL0(dir >= 0 && dir <= 2, "input dir is out of range");
int i;
int order0 = m_base[0]->GetNumModes ();
int order1 = m_base[1]->GetNumModes ();
int nquad0 = m_base[0]->GetNumPoints();
int nquad1 = m_base[1]->GetNumPoints();
int nquad2 = m_base[2]->GetNumPoints();
const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
const Array<OneD, const NekDouble> &z2 = m_base[2]->GetZ();
Array<OneD, NekDouble> gfac0(nquad0);
Array<OneD, NekDouble> gfac2(nquad2);
Array<OneD, NekDouble> tmp0 (nquad0*nquad1*nquad2);
Array<OneD, NekDouble> wsp (order0*nquad2*(nquad1+order1));
// set up geometric factor: (1+z0)/2
for (i = 0; i < nquad0; ++i)
{
gfac0[i] = 0.5*(1+z0[i]);
}
// Set up geometric factor: 2/(1-z2)
for (i = 0; i < nquad2; ++i)
{
gfac2[i] = 2.0/(1-z2[i]);
}
// Scale first derivative term by gfac2.
if (dir != 1)
{
for (i = 0; i < nquad2; ++i)
{
Vmath::Smul(nquad0*nquad1,gfac2[i],
&inarray[0]+i*nquad0*nquad1,1,
&tmp0 [0]+i*nquad0*nquad1,1);
}
}
switch (dir)
{
case 0:
{
m_base[1]->GetBdata (),
m_base[2]->GetBdata (),
tmp0,outarray,wsp,
true,true,true);
break;
}
case 1:
{
m_base[1]->GetDbdata(),
m_base[2]->GetBdata (),
tmp0,outarray,wsp,
true,true,true);
break;
}
case 2:
{
Array<OneD, NekDouble> tmp1(m_ncoeffs);
// Scale eta_1 derivative with gfac0.
for(i = 0; i < nquad1*nquad2; ++i)
{
Vmath::Vmul(nquad0,&gfac0[0],1,&tmp0[0]+i*nquad0,1,&tmp0[0]+i*nquad0,1);
}
m_base[1]->GetBdata(),
m_base[2]->GetBdata(),
tmp0,tmp1,wsp,
true,true,true);
m_base[1]->GetBdata(),
m_base[2]->GetDbdata(),
tmp0,outarray,wsp,
true,true,true);
Vmath::Vadd(m_ncoeffs,&tmp1[0],1,&outarray[0],1,&outarray[0],1);
break;
}
}
}
bool Nektar::StdRegions::StdPrismExp::v_IsBoundaryInteriorExpansion ( )
protectedvirtual
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 664 of file StdPrismExp.cpp.

References Nektar::NekConstants::kNekZeroTol.

{
if( fabs(xi[2]-1.0) < NekConstants::kNekZeroTol)
{
// Very top point of the prism
eta[0] = -1.0;
eta[1] = xi[1];
eta[2] = 1.0;
}
else
{
// Third basis function collapsed to "pr" direction instead of
// "qr" direction
eta[2] = xi[2]; // eta_z = xi_z
eta[1] = xi[1]; //eta_y = xi_y
eta[0] = 2.0*(1.0 + xi[0])/(1.0 - xi[2]) - 1.0;
}
}
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 1811 of file StdPrismExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 766 of file StdPrismExp.cpp.

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

{
"BasisType is not a boundary interior form");
"BasisType is not a boundary interior form");
"BasisType is not a boundary interior form");
int P = m_base[0]->GetNumModes()-1;
int Q = m_base[1]->GetNumModes()-1;
int R = m_base[2]->GetNumModes()-1;
return (P+1)*(Q+1) // 1 rect. face on base
+ 2*(Q+1)*(R+1) // other 2 rect. faces
+ 2*(R+1) + P*(1 + 2*R - P); // 2 tri. faces
}
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 96 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().

{
int Qx = m_base[0]->GetNumPoints();
int Qy = m_base[1]->GetNumPoints();
int Qz = m_base[2]->GetNumPoints();
int Qtot = Qx*Qy*Qz;
Array<OneD, NekDouble> dEta_bar1(Qtot,0.0);
Array<OneD, const NekDouble> eta_x, eta_z;
eta_x = m_base[0]->GetZ();
eta_z = m_base[2]->GetZ();
int i, k;
bool Do_1 = (out_dxi1.num_elements() > 0)? true:false;
bool Do_3 = (out_dxi3.num_elements() > 0)? true:false;
// out_dXi2 is just a tensor derivative so is just passed through
if(Do_3)
{
PhysTensorDeriv(u_physical, dEta_bar1, out_dxi2, out_dxi3);
}
else if(Do_1)
{
PhysTensorDeriv(u_physical, dEta_bar1, out_dxi2, NullNekDouble1DArray);
}
else // case if just require 2nd direction
{
}
if(Do_1)
{
for (k = 0; k < Qz; ++k)
{
Vmath::Smul(Qx*Qy,2.0/(1.0-eta_z[k]),&dEta_bar1[0] + k*Qx*Qy,1,
&out_dxi1[0] + k*Qx*Qy,1);
}
}
if(Do_3)
{
// divide dEta_Bar1 by (1-eta_z)
for (k = 0; k < Qz; ++k)
{
Vmath::Smul(Qx*Qy, 1.0/(1.0-eta_z[k]),&dEta_bar1[0]+k*Qx*Qy,1,
&dEta_bar1[0]+k*Qx*Qy,1);
}
// Multiply dEta_Bar1 by (1+eta_x) and add ot out_dxi3
for (i = 0; i < Qx; ++i)
{
Vmath::Svtvp(Qz*Qy,1.0+eta_x[i],&dEta_bar1[0]+i,Qx,
&out_dxi3[0]+i,Qx,&out_dxi3[0]+i,Qx);
}
}
}
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 160 of file StdPrismExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1931 of file StdPrismExp.cpp.

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

{
int nquad0 = m_base[0]->GetNumPoints();
int nquad1 = m_base[1]->GetNumPoints();
int nquad2 = m_base[2]->GetNumPoints();
int nqtot = nquad0*nquad1*nquad2;
int nmodes0 = m_base[0]->GetNumModes();
int nmodes1 = m_base[1]->GetNumModes();
int nmodes2 = m_base[2]->GetNumModes();
int numMax = nmodes0;
Array<OneD, NekDouble> coeff (m_ncoeffs);
Array<OneD, NekDouble> coeff_tmp1(m_ncoeffs, 0.0);
Array<OneD, NekDouble> phys_tmp (nqtot, 0.0);
Array<OneD, NekDouble> tmp, tmp2, tmp3, tmp4;
const LibUtilities::PointsKey Pkey0 = m_base[0]->GetPointsKey();
const LibUtilities::PointsKey Pkey1 = m_base[1]->GetPointsKey();
const LibUtilities::PointsKey Pkey2 = m_base[2]->GetPointsKey();
LibUtilities::BasisKey bortho0(
LibUtilities::eOrtho_A, nmodes0, Pkey0);
LibUtilities::BasisKey bortho1(
LibUtilities::eOrtho_A, nmodes1, Pkey1);
LibUtilities::BasisKey bortho2(
LibUtilities::eOrtho_B, nmodes2, Pkey2);
int cnt = 0;
int u = 0;
int i = 0;
OrthoPrismExp = MemoryManager<StdRegions::StdPrismExp>
::AllocateSharedPtr(bortho0, bortho1, bortho2);
BwdTrans(inarray,phys_tmp);
OrthoPrismExp->FwdTrans(phys_tmp, coeff);
// filtering
for (u = 0; u < numMin; ++u)
{
for (i = 0; i < numMin; ++i)
{
Vmath::Vcopy(numMin - u, tmp = coeff + cnt, 1,
tmp2 = coeff_tmp1 + cnt, 1);
cnt += numMax - u;
}
for (i = numMin; i < numMax; ++i)
{
cnt += numMax - u;
}
}
OrthoPrismExp->BwdTrans(coeff_tmp1, phys_tmp);
StdPrismExp::FwdTrans(phys_tmp, outarray);
}
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 195 of file StdPrismExp.cpp.

References v_PhysDeriv().

{
StdPrismExp::v_PhysDeriv(inarray, out_d0, out_d1, out_d2);
}
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 203 of file StdPrismExp.cpp.

References v_PhysDeriv().

{
StdPrismExp::v_PhysDeriv(dir, inarray, outarray);
}
void Nektar::StdRegions::StdPrismExp::v_SVVLaplacianFilter ( Array< OneD, NekDouble > &  array,
const StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

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

{
// Generate an orthonogal expansion
int qa = m_base[0]->GetNumPoints();
int qb = m_base[1]->GetNumPoints();
int qc = m_base[2]->GetNumPoints();
int nmodes_a = m_base[0]->GetNumModes();
int nmodes_b = m_base[1]->GetNumModes();
int nmodes_c = m_base[2]->GetNumModes();
// Declare orthogonal basis.
LibUtilities::PointsKey pa(qa,m_base[0]->GetPointsType());
LibUtilities::PointsKey pb(qb,m_base[1]->GetPointsType());
LibUtilities::PointsKey pc(qc,m_base[2]->GetPointsType());
LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A,nmodes_a,pa);
LibUtilities::BasisKey Bb(LibUtilities::eOrtho_A,nmodes_b,pb);
LibUtilities::BasisKey Bc(LibUtilities::eOrtho_B,nmodes_c,pc);
StdPrismExp OrthoExp(Ba,Bb,Bc);
Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
int i,j,k,cnt = 0;
//SVV filter paramaters (how much added diffusion relative to physical one
// and fraction of modes from which you start applying this added diffusion)
//
NekDouble SvvDiffCoeff = mkey.GetConstFactor(StdRegions::eFactorSVVDiffCoeff);
NekDouble SVVCutOff = mkey.GetConstFactor(StdRegions::eFactorSVVCutoffRatio);
//Defining the cut of mode
int cutoff_a = (int) (SVVCutOff*nmodes_a);
int cutoff_b = (int) (SVVCutOff*nmodes_b);
int cutoff_c = (int) (SVVCutOff*nmodes_c);
//To avoid the fac[j] from blowing up
NekDouble epsilon = 1;
// project onto modal space.
OrthoExp.FwdTrans(array,orthocoeffs);
int nmodes = min(min(nmodes_a,nmodes_b),nmodes_c);
NekDouble cutoff = min(min(cutoff_a,cutoff_b),cutoff_c);
//------"New" Version August 22nd '13--------------------
for(i = 0; i < nmodes_a; ++i)//P
{
for(j = 0; j < nmodes_b; ++j) //Q
{
for(k = 0; k < nmodes_c-i; ++k) //R
{
if(j >= cutoff || i + k >= cutoff)
{
orthocoeffs[cnt] *= (1.0+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)))));
}
cnt++;
}
}
}
// backward transform to physical space
OrthoExp.BwdTrans(orthocoeffs,array);
}