Nektar++
Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes | List of all members
Nektar::LocalRegions::PyrExp Class Reference

#include <PyrExp.h>

Inheritance diagram for Nektar::LocalRegions::PyrExp:
Inheritance graph
[legend]
Collaboration diagram for Nektar::LocalRegions::PyrExp:
Collaboration graph
[legend]

Public Member Functions

 PyrExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc, const SpatialDomains::PyrGeomSharedPtr &geom)
 Constructor using BasisKey class for quadrature points and order definition. More...
 
 PyrExp (const PyrExp &T)
 
 ~PyrExp ()
 
- Public Member Functions inherited from Nektar::StdRegions::StdPyrExp
 StdPyrExp ()
 
 StdPyrExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
 
 StdPyrExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc, NekDouble *coeffs, NekDouble *phys)
 
 StdPyrExp (const StdPyrExp &T)
 
 ~StdPyrExp ()
 
int GetTetMode (int I, int J, int K)
 Number tetrahedral modes in r-direction. Much the same as StdTetExp::GetTetMode but slightly simplified since we know that the polynomial order is the same in each direction. More...
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion3D
 StdExpansion3D ()
 
 StdExpansion3D (int numcoeffs, const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
 
 StdExpansion3D (const StdExpansion3D &T)
 
virtual ~StdExpansion3D ()
 
void PhysTensorDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d1, Array< OneD, NekDouble > &outarray_d2, Array< OneD, NekDouble > &outarray_d3)
 Calculate the 3D derivative in the local tensor/collapsed coordinate at the physical points. More...
 
void BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
 
void IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion
 StdExpansion ()
 Default Constructor. More...
 
 StdExpansion (const int numcoeffs, const int numbases, const LibUtilities::BasisKey &Ba=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bb=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bc=LibUtilities::NullBasisKey)
 Constructor. More...
 
 StdExpansion (const StdExpansion &T)
 Copy Constructor. More...
 
virtual ~StdExpansion ()
 Destructor. More...
 
int GetNumBases () const
 This function returns the number of 1D bases used in the expansion. More...
 
const Array< OneD, const LibUtilities::BasisSharedPtr > & GetBase () const
 This function gets the shared point to basis. More...
 
const LibUtilities::BasisSharedPtrGetBasis (int dir) const
 This function gets the shared point to basis in the dir direction. More...
 
int GetNcoeffs (void) const
 This function returns the total number of coefficients used in the expansion. More...
 
int GetTotPoints () const
 This function returns the total number of quadrature points used in the element. More...
 
LibUtilities::BasisType GetBasisType (const int dir) const
 This function returns the type of basis used in the dir direction. More...
 
int GetBasisNumModes (const int dir) const
 This function returns the number of expansion modes in the dir direction. More...
 
int EvalBasisNumModesMax (void) const
 This function returns the maximum number of expansion modes over all local directions. More...
 
LibUtilities::PointsType GetPointsType (const int dir) const
 This function returns the type of quadrature points used in the dir direction. More...
 
int GetNumPoints (const int dir) const
 This function returns the number of quadrature points in the dir direction. More...
 
const Array< OneD, const NekDouble > & GetPoints (const int dir) const
 This function returns a pointer to the array containing the quadrature points in dir direction. More...
 
int GetNverts () const
 This function returns the number of vertices of the expansion domain. More...
 
int 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
 
LibUtilities::PointsKey GetFacePointsKey (const int i, const int j) const
 
int NumBndryCoeffs (void) const
 
int NumDGBndryCoeffs (void) const
 
LibUtilities::BasisType GetEdgeBasisType (const int i) const
 This function returns the type of expansion basis on the i-th edge. More...
 
const LibUtilities::PointsKey GetNodalPointsKey () const
 This function returns the type of expansion Nodal point type if defined. More...
 
int GetNfaces () const
 This function returns the number of faces of the expansion domain. More...
 
int GetNtrace () const
 Returns the number of trace elements connected to this element. More...
 
LibUtilities::ShapeType DetShapeType () const
 This function returns the shape of the expansion domain. More...
 
boost::shared_ptr< StdExpansionGetStdExp (void) const
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion ()
 
bool IsNodalNonTensorialExp ()
 
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space. More...
 
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space. More...
 
void FwdTrans_BndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble Integral (const Array< OneD, const NekDouble > &inarray)
 This function integrates the specified function over the domain. More...
 
void FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 This function fills the array outarray with the mode-th mode of the expansion. More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 this function calculates the inner product of a given function f with the different modes of the expansion More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check)
 
void IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
int GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void SetElmtId (const int id)
 Set the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2=NullNekDouble1DArray, Array< OneD, NekDouble > &coords_3=NullNekDouble1DArray)
 this function returns the physical coordinates of the quadrature points of the expansion More...
 
void GetCoord (const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
 given the coordinates of a point of the element in the local collapsed coordinate system, this function calculates the physical coordinates of the point More...
 
DNekMatSharedPtr GetStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr GetStdStaticCondMatrix (const StdMatrixKey &mkey)
 
IndexMapValuesSharedPtr GetIndexMap (const IndexMapKey &ikey)
 
const Array< OneD, const NekDouble > & GetPhysNormals (void)
 
void SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
virtual void SetUpPhysNormals (const int edge)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, 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. 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 MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
DNekMatSharedPtr CreateGeneralMatrix (const StdMatrixKey &mkey)
 this function generates the mass matrix $\mathbf{M}[i][j] = \int \phi_i(\mathbf{x}) \phi_j(\mathbf{x}) d\mathbf{x}$ More...
 
void GeneralMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
DNekMatSharedPtr GenMatrix (const StdMatrixKey &mkey)
 
void PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void PhysDeriv_s (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds)
 
void PhysDeriv_n (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn)
 
void PhysDirectionalDeriv (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &outarray)
 
void StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void AddRobinMassMatrix (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
void AddRobinEdgeContribution (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, Array< OneD, NekDouble > &coeffs)
 
NekDouble PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta. More...
 
const boost::shared_ptr< SpatialDomains::GeomFactors > & GetMetricInfo (void) const
 
virtual int v_GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
virtual const Array< OneD, const NekDouble > & v_GetPhysNormals (void)
 
virtual void v_SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
virtual void v_SetUpPhysNormals (const int edge)
 
virtual void v_ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmode_offset, NekDouble *coeffs)
 Unpack data from input file assuming it comes from the same expansion type. More...
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, 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 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)
 
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. 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...
 
template<class T >
boost::shared_ptr< T > as ()
 
void IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
- Public Member Functions inherited from Nektar::LocalRegions::Expansion3D
 Expansion3D (SpatialDomains::Geometry3DSharedPtr pGeom)
 
virtual ~Expansion3D ()
 
void SetFaceExp (const int face, Expansion2DSharedPtr &f)
 
Expansion2DSharedPtr GetFaceExp (const int face)
 
void SetTraceToGeomOrientation (Array< OneD, NekDouble > &inout)
 Align trace orientation with the geometry orientation. More...
 
void SetFaceToGeomOrientation (const int face, Array< OneD, NekDouble > &inout)
 Align face orientation with the geometry orientation. More...
 
void AddHDGHelmholtzFaceTerms (const NekDouble tau, const int edge, Array< OneD, NekDouble > &facePhys, const StdRegions::VarCoeffMap &dirForcing, Array< OneD, NekDouble > &outarray)
 
void AddNormTraceInt (const int dir, Array< OneD, ExpansionSharedPtr > &FaceExp, Array< OneD, Array< OneD, NekDouble > > &faceCoeffs, Array< OneD, NekDouble > &outarray)
 
void AddNormTraceInt (const int dir, Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &FaceExp, Array< OneD, NekDouble > &outarray, const StdRegions::VarCoeffMap &varcoeffs)
 
void AddFaceBoundaryInt (const int face, ExpansionSharedPtr &FaceExp, Array< OneD, NekDouble > &facePhys, Array< OneD, NekDouble > &outarray, const StdRegions::VarCoeffMap &varcoeffs=StdRegions::NullVarCoeffMap)
 
SpatialDomains::Geometry3DSharedPtr GetGeom3D () const
 
- Public Member Functions inherited from Nektar::LocalRegions::Expansion
 Expansion (SpatialDomains::GeometrySharedPtr pGeom)
 
 Expansion (const Expansion &pSrc)
 
virtual ~Expansion ()
 
DNekScalMatSharedPtr GetLocMatrix (const LocalRegions::MatrixKey &mkey)
 
DNekScalMatSharedPtr GetLocMatrix (const StdRegions::MatrixType mtype, const StdRegions::ConstFactorMap &factors=StdRegions::NullConstFactorMap, const StdRegions::VarCoeffMap &varcoeffs=StdRegions::NullVarCoeffMap)
 
SpatialDomains::GeometrySharedPtr GetGeom () const
 
void Reset ()
 
virtual const SpatialDomains::GeomFactorsSharedPtrv_GetMetricInfo () const
 
DNekMatSharedPtr BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
 
DNekMatSharedPtr BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd)
 
void AddEdgeNormBoundaryInt (const int edge, const boost::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
void AddEdgeNormBoundaryInt (const int edge, const boost::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
void AddFaceNormBoundaryInt (const int face, const boost::shared_ptr< Expansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
void DGDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &coeffs, Array< OneD, NekDouble > &outarray)
 

Protected Member Functions

virtual NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray)
 Integrate the physical point list inarray over pyramidic region and return the value. More...
 
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_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 (this)->m_coeffs. 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 StdRegions::StdExpansionSharedPtr v_GetStdExp (void) const
 
virtual void v_GetCoord (const Array< OneD, const NekDouble > &Lcoords, Array< OneD, NekDouble > &coords)
 
virtual void v_GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3)
 
virtual NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &coord, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
virtual int v_GetCoordim ()
 
void v_ComputeFaceNormal (const int face)
 
virtual void v_GetFacePhysVals (const int face, const StdRegions::StdExpansionSharedPtr &FaceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient)
 
virtual DNekMatSharedPtr v_GenMatrix (const StdRegions::StdMatrixKey &mkey)
 
virtual DNekMatSharedPtr v_CreateStdMatrix (const StdRegions::StdMatrixKey &mkey)
 
virtual DNekScalMatSharedPtr v_GetLocMatrix (const MatrixKey &mkey)
 
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const MatrixKey &mkey)
 
void v_DropLocStaticCondMatrix (const MatrixKey &mkey)
 
DNekScalMatSharedPtr CreateMatrix (const MatrixKey &mkey)
 
DNekScalBlkMatSharedPtr CreateStaticCondMatrix (const MatrixKey &mkey)
 
virtual void v_ComputeLaplacianMetric ()
 
- Protected Member Functions inherited from Nektar::StdRegions::StdPyrExp
virtual void v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Calculate the derivative of the physical points in a given direction. More...
 
virtual void v_StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
 
virtual void v_StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Backward transformation is evaluated at the quadrature points. More...
 
virtual void v_BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
 
virtual void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
virtual void v_IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
 
virtual void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 
virtual void v_FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 
virtual int v_GetNverts () const
 
virtual int v_GetNedges () const
 
virtual int v_GetNfaces () const
 
virtual LibUtilities::ShapeType v_DetShapeType () const
 
virtual int v_NumBndryCoeffs () const
 
virtual int v_GetEdgeNcoeffs (const int i) const
 
virtual int v_GetFaceNcoeffs (const int i) const
 
virtual int v_GetFaceIntNcoeffs (const int i) const
 
virtual int v_GetFaceNumPoints (const int i) const
 
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
virtual const LibUtilities::BasisKey v_DetFaceBasisKey (const int i, const int k) const
 
virtual LibUtilities::BasisType v_GetEdgeBasisType (const int i) const
 
virtual void v_GetFaceToElementMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int nummodesA=-1, int nummodesB=-1)
 
virtual int v_GetVertexMap (int localVertexId, bool useCoeffPacking=false)
 
virtual void v_GetEdgeInteriorMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
 
virtual void v_GetFaceInteriorMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
 
virtual void v_GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
virtual void v_GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion3D
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 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. 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)
 
- Protected Member Functions inherited from Nektar::LocalRegions::Expansion3D
virtual void v_DGDeriv (const int dir, const Array< OneD, const NekDouble > &incoeffs, Array< OneD, ExpansionSharedPtr > &FaceExp, Array< OneD, Array< OneD, NekDouble > > &faceCoeffs, Array< OneD, NekDouble > &out_d)
 Evaluate coefficients of weak deriviative in the direction dir given the input coefficicents incoeffs and the imposed boundary values in EdgeExp (which will have its phys space updated). More...
 
virtual void v_AddFaceNormBoundaryInt (const int face, const ExpansionSharedPtr &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddRobinMassMatrix (const int face, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
virtual StdRegions::Orientation v_GetForient (int face)
 
virtual Array< OneD, unsigned int > v_GetEdgeInverseBoundaryMap (int eid)
 
virtual Array< OneD, unsigned int > v_GetFaceInverseBoundaryMap (int fid, StdRegions::Orientation faceOrient=StdRegions::eNoOrientation)
 
virtual DNekMatSharedPtr v_BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
 
virtual DNekMatSharedPtr v_BuildInverseTransformationMatrix (const DNekScalMatSharedPtr &m_transformationmatrix)
 Build inverse and inverse transposed transformation matrix: $\mathbf{R^{-1}}$ and $\mathbf{R^{-T}}$. More...
 
virtual DNekMatSharedPtr v_BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd)
 
- Protected Member Functions inherited from Nektar::LocalRegions::Expansion
void ComputeLaplacianMetric ()
 
void ComputeQuadratureMetric ()
 
virtual void v_MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddEdgeNormBoundaryInt (const int edge, const boost::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddEdgeNormBoundaryInt (const int edge, const boost::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddFaceNormBoundaryInt (const int face, const boost::shared_ptr< Expansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 

Private Member Functions

virtual void v_LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
 

Private Attributes

LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLessm_matrixManager
 
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLessm_staticCondMatrixManager
 

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::BasisSharedPtrm_base
 
int m_elmt_id
 
int m_ncoeffs
 
LibUtilities::NekManager< StdMatrixKey, DNekMat, StdMatrixKey::opLessm_stdMatrixManager
 
LibUtilities::NekManager< StdMatrixKey, DNekBlkMat, StdMatrixKey::opLessm_stdStaticCondMatrixManager
 
LibUtilities::NekManager< IndexMapKey, IndexMapValues, IndexMapKey::opLessm_IndexMapManager
 
- Protected Attributes inherited from Nektar::LocalRegions::Expansion
SpatialDomains::GeometrySharedPtr m_geom
 
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
 
MetricMap m_metrics
 

Detailed Description

Definition at line 51 of file PyrExp.h.

Constructor & Destructor Documentation

Nektar::LocalRegions::PyrExp::PyrExp ( const LibUtilities::BasisKey Ba,
const LibUtilities::BasisKey Bb,
const LibUtilities::BasisKey Bc,
const SpatialDomains::PyrGeomSharedPtr geom 
)

Constructor using BasisKey class for quadrature points and order definition.

Definition at line 44 of file PyrExp.cpp.

47  :
49  Ba.GetNumModes(),
50  Bb.GetNumModes(),
51  Bc.GetNumModes()),
52  3, Ba, Bb, Bc),
54  Ba.GetNumModes(),
55  Bb.GetNumModes(),
56  Bc.GetNumModes()),
57  Ba, Bb, Bc),
58  StdPyrExp (Ba,Bb,Bc),
59  Expansion (geom),
60  Expansion3D (geom),
62  boost::bind(&PyrExp::CreateMatrix, this, _1),
63  std::string("PyrExpMatrix")),
65  boost::bind(&PyrExp::CreateStaticCondMatrix, this, _1),
66  std::string("PyrExpStaticCondMatrix"))
67  {
68  }
DNekScalBlkMatSharedPtr CreateStaticCondMatrix(const MatrixKey &mkey)
Definition: PyrExp.cpp:957
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: PyrExp.h:154
DNekScalMatSharedPtr CreateMatrix(const MatrixKey &mkey)
Definition: PyrExp.cpp:831
Expansion(SpatialDomains::GeometrySharedPtr pGeom)
Definition: Expansion.cpp:46
int getNumberOfCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:232
StdExpansion()
Default Constructor.
Expansion3D(SpatialDomains::Geometry3DSharedPtr pGeom)
Definition: Expansion3D.h:63
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: PyrExp.h:153
Nektar::LocalRegions::PyrExp::PyrExp ( const PyrExp T)

Definition at line 70 of file PyrExp.cpp.

70  :
71  StdExpansion (T),
72  StdExpansion3D(T),
73  StdPyrExp (T),
74  Expansion (T),
75  Expansion3D (T),
76  m_matrixManager(T.m_matrixManager),
77  m_staticCondMatrixManager(T.m_staticCondMatrixManager)
78  {
79  }
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: PyrExp.h:154
Expansion(SpatialDomains::GeometrySharedPtr pGeom)
Definition: Expansion.cpp:46
StdExpansion()
Default Constructor.
Expansion3D(SpatialDomains::Geometry3DSharedPtr pGeom)
Definition: Expansion3D.h:63
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: PyrExp.h:153
Nektar::LocalRegions::PyrExp::~PyrExp ( )

Definition at line 81 of file PyrExp.cpp.

82  {
83  }

Member Function Documentation

DNekScalMatSharedPtr Nektar::LocalRegions::PyrExp::CreateMatrix ( const MatrixKey mkey)
protected

Definition at line 831 of file PyrExp.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL2, Nektar::StdRegions::StdMatrixKey::ConstFactorExists(), Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eFactorSVVCutoffRatio, ErrorUtil::efatal, Nektar::StdRegions::eHelmholtz, Nektar::StdRegions::eInvMass, Nektar::StdRegions::eLaplacian, Nektar::StdRegions::eLaplacian00, Nektar::StdRegions::eLaplacian01, Nektar::StdRegions::eLaplacian02, Nektar::StdRegions::eLaplacian11, Nektar::StdRegions::eLaplacian12, Nektar::StdRegions::eLaplacian22, Nektar::StdRegions::eMass, Nektar::SpatialDomains::eNoGeomType, Nektar::StdRegions::StdExpansion::GenMatrix(), Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdMatrixKey::GetConstFactors(), Nektar::StdRegions::StdMatrixKey::GetMatrixType(), Nektar::StdRegions::StdMatrixKey::GetNVarCoeff(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdMatrixKey::GetShapeType(), Nektar::StdRegions::StdExpansion::GetStdMatrix(), Nektar::StdRegions::StdMatrixKey::GetVarCoeffs(), m_matrixManager, Nektar::LocalRegions::Expansion::m_metricinfo, NEKERROR, and Nektar::Transpose().

832  {
833  DNekScalMatSharedPtr returnval;
835 
836  ASSERTL2(m_metricinfo->GetGtype() != SpatialDomains::eNoGeomType,"Geometric information is not set up");
837 
838  switch(mkey.GetMatrixType())
839  {
840  case StdRegions::eMass:
841  {
842  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
843  {
844  NekDouble one = 1.0;
845  DNekMatSharedPtr mat = GenMatrix(mkey);
847  }
848  else
849  {
850  NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
851  DNekMatSharedPtr mat = GetStdMatrix(mkey);
853  }
854  }
855  break;
857  {
858  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
859  {
860  NekDouble one = 1.0;
861  StdRegions::StdMatrixKey masskey(StdRegions::eMass,DetShapeType(),
862  *this);
863  DNekMatSharedPtr mat = GenMatrix(masskey);
864  mat->Invert();
866  }
867  else
868  {
869  NekDouble fac = 1.0/(m_metricinfo->GetJac(ptsKeys))[0];
870  DNekMatSharedPtr mat = GetStdMatrix(mkey);
872  }
873  }
874  break;
876  {
877  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed ||
878  mkey.GetNVarCoeff() > 0 ||
879  mkey.ConstFactorExists(
881  {
882  NekDouble one = 1.0;
883  DNekMatSharedPtr mat = GenMatrix(mkey);
884 
886  }
887  else
888  {
889  MatrixKey lap00key(StdRegions::eLaplacian00,
890  mkey.GetShapeType(), *this);
891  MatrixKey lap01key(StdRegions::eLaplacian01,
892  mkey.GetShapeType(), *this);
893  MatrixKey lap02key(StdRegions::eLaplacian02,
894  mkey.GetShapeType(), *this);
895  MatrixKey lap11key(StdRegions::eLaplacian11,
896  mkey.GetShapeType(), *this);
897  MatrixKey lap12key(StdRegions::eLaplacian12,
898  mkey.GetShapeType(), *this);
899  MatrixKey lap22key(StdRegions::eLaplacian22,
900  mkey.GetShapeType(), *this);
901 
902  DNekMat &lap00 = *GetStdMatrix(lap00key);
903  DNekMat &lap01 = *GetStdMatrix(lap01key);
904  DNekMat &lap02 = *GetStdMatrix(lap02key);
905  DNekMat &lap11 = *GetStdMatrix(lap11key);
906  DNekMat &lap12 = *GetStdMatrix(lap12key);
907  DNekMat &lap22 = *GetStdMatrix(lap22key);
908 
909  NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
910  Array<TwoD, const NekDouble> gmat =
911  m_metricinfo->GetGmat(ptsKeys);
912 
913  int rows = lap00.GetRows();
914  int cols = lap00.GetColumns();
915 
917  ::AllocateSharedPtr(rows,cols);
918 
919  (*lap) = gmat[0][0]*lap00
920  + gmat[4][0]*lap11
921  + gmat[8][0]*lap22
922  + gmat[3][0]*(lap01 + Transpose(lap01))
923  + gmat[6][0]*(lap02 + Transpose(lap02))
924  + gmat[7][0]*(lap12 + Transpose(lap12));
925 
926  returnval = MemoryManager<DNekScalMat>
927  ::AllocateSharedPtr(jac, lap);
928  }
929  }
930  break;
932  {
933  NekDouble factor = mkey.GetConstFactor(StdRegions::eFactorLambda);
934  MatrixKey masskey(StdRegions::eMass, mkey.GetShapeType(), *this);
935  DNekScalMat &MassMat = *(this->m_matrixManager[masskey]);
936  MatrixKey lapkey(StdRegions::eLaplacian, mkey.GetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
937  DNekScalMat &LapMat = *(this->m_matrixManager[lapkey]);
938 
939  int rows = LapMat.GetRows();
940  int cols = LapMat.GetColumns();
941 
943 
944  (*helm) = LapMat + factor*MassMat;
945 
946  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(1.0, helm);
947  }
948  break;
949  default:
950  NEKERROR(ErrorUtil::efatal, "Matrix creation not defined");
951  break;
952  }
953 
954  return returnval;
955  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:454
DNekMatSharedPtr GenMatrix(const StdMatrixKey &mkey)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:158
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:220
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:126
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:684
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
NekMatrix< InnerMatrixType, BlockMatrixTag > Transpose(NekMatrix< InnerMatrixType, BlockMatrixTag > &rhs)
NekMatrix< NekDouble, StandardMatrixTag > DNekMat
Definition: NekTypeDefs.hpp:52
double NekDouble
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:213
Geometry is curved or has non-constant factors.
NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag > DNekScalMat
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: PyrExp.h:153
DNekScalBlkMatSharedPtr Nektar::LocalRegions::PyrExp::CreateStaticCondMatrix ( const MatrixKey mkey)
protected

Definition at line 957 of file PyrExp.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL2, Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::eHelmholtz, Nektar::StdRegions::eLaplacian, Nektar::SpatialDomains::eNoGeomType, Nektar::StdRegions::StdExpansion::GetBoundaryMap(), Nektar::StdRegions::StdExpansion::GetInteriorMap(), Nektar::LocalRegions::Expansion::GetLocMatrix(), Nektar::StdRegions::StdMatrixKey::GetMatrixType(), Nektar::StdRegions::StdExpansion::GetStdStaticCondMatrix(), Nektar::LocalRegions::Expansion::m_metricinfo, Nektar::StdRegions::StdExpansion::m_ncoeffs, and Nektar::StdRegions::StdExpansion::NumBndryCoeffs().

958  {
959  DNekScalBlkMatSharedPtr returnval;
960 
961  ASSERTL2(m_metricinfo->GetGtype() != SpatialDomains::eNoGeomType,"Geometric information is not set up");
962 
963  // set up block matrix system
964  unsigned int nbdry = NumBndryCoeffs();
965  unsigned int nint = (unsigned int)(m_ncoeffs - nbdry);
966  unsigned int exp_size[] = {nbdry, nint};
967  unsigned int nblks = 2;
968  returnval = MemoryManager<DNekScalBlkMat>::AllocateSharedPtr(nblks, nblks, exp_size, exp_size); //Really need a constructor which takes Arrays
969  NekDouble factor = 1.0;
970 
971  switch(mkey.GetMatrixType())
972  {
974  case StdRegions::eHelmholtz: // special case since Helmholtz not defined in StdRegions
975 
976  // use Deformed case for both regular and deformed geometries
977  factor = 1.0;
978  goto UseLocRegionsMatrix;
979  break;
980  default:
981  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
982  {
983  factor = 1.0;
984  goto UseLocRegionsMatrix;
985  }
986  else
987  {
989  factor = mat->Scale();
990  goto UseStdRegionsMatrix;
991  }
992  break;
993  UseStdRegionsMatrix:
994  {
995  NekDouble invfactor = 1.0/factor;
996  NekDouble one = 1.0;
999  DNekMatSharedPtr Asubmat;
1000 
1001  //TODO: check below
1002  returnval->SetBlock(0,0,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(factor,Asubmat = mat->GetBlock(0,0)));
1003  returnval->SetBlock(0,1,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,Asubmat = mat->GetBlock(0,1)));
1004  returnval->SetBlock(1,0,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(factor,Asubmat = mat->GetBlock(1,0)));
1005  returnval->SetBlock(1,1,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(invfactor,Asubmat = mat->GetBlock(1,1)));
1006  }
1007  break;
1008  UseLocRegionsMatrix:
1009  {
1010  int i,j;
1011  NekDouble invfactor = 1.0/factor;
1012  NekDouble one = 1.0;
1013  DNekScalMat &mat = *GetLocMatrix(mkey);
1018 
1019  Array<OneD,unsigned int> bmap(nbdry);
1020  Array<OneD,unsigned int> imap(nint);
1021  GetBoundaryMap(bmap);
1022  GetInteriorMap(imap);
1023 
1024  for(i = 0; i < nbdry; ++i)
1025  {
1026  for(j = 0; j < nbdry; ++j)
1027  {
1028  (*A)(i,j) = mat(bmap[i],bmap[j]);
1029  }
1030 
1031  for(j = 0; j < nint; ++j)
1032  {
1033  (*B)(i,j) = mat(bmap[i],imap[j]);
1034  }
1035  }
1036 
1037  for(i = 0; i < nint; ++i)
1038  {
1039  for(j = 0; j < nbdry; ++j)
1040  {
1041  (*C)(i,j) = mat(imap[i],bmap[j]);
1042  }
1043 
1044  for(j = 0; j < nint; ++j)
1045  {
1046  (*D)(i,j) = mat(imap[i],imap[j]);
1047  }
1048  }
1049 
1050  // Calculate static condensed system
1051  if(nint)
1052  {
1053  D->Invert();
1054  (*B) = (*B)*(*D);
1055  (*A) = (*A) - (*B)*(*C);
1056  }
1057 
1058  DNekScalMatSharedPtr Atmp;
1059 
1060  returnval->SetBlock(0,0,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(factor,A));
1061  returnval->SetBlock(0,1,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,B));
1062  returnval->SetBlock(1,0,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(factor,C));
1063  returnval->SetBlock(1,1,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(invfactor,D));
1064 
1065  }
1066  }
1067  return returnval;
1068  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:126
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
DNekBlkMatSharedPtr GetStdStaticCondMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:689
boost::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:74
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:795
double NekDouble
boost::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:72
DNekScalMatSharedPtr GetLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:83
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:213
Geometry is curved or has non-constant factors.
NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag > DNekScalMat
void GetBoundaryMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:790
void Nektar::LocalRegions::PyrExp::v_ComputeFaceNormal ( const int  face)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 512 of file PyrExp.cpp.

References ASSERTL0, Nektar::StdRegions::StdExpansion::DetFaceBasisKey(), Nektar::SpatialDomains::eMovingRegular, Nektar::SpatialDomains::eRegular, Vmath::Fill(), Nektar::StdRegions::StdExpansion::GetCoordim(), Nektar::LocalRegions::Expansion::GetGeom(), Nektar::LibUtilities::BasisKey::GetNumPoints(), Nektar::LibUtilities::BasisKey::GetPointsKey(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::LibUtilities::Interp2D(), Nektar::StdRegions::StdExpansion3D::m_faceNormals, Vmath::Sdiv(), Vmath::Vmul(), Vmath::Vsqrt(), Vmath::Vvtvp(), and Vmath::Zero().

513  {
514  const SpatialDomains::GeomFactorsSharedPtr &geomFactors =
515  GetGeom()->GetMetricInfo();
517  SpatialDomains::GeomType type = geomFactors->GetGtype();
518  const Array<TwoD, const NekDouble> &df = geomFactors->GetDerivFactors(ptsKeys);
519  const Array<OneD, const NekDouble> &jac = geomFactors->GetJac(ptsKeys);
520 
521  LibUtilities::BasisKey tobasis0 = DetFaceBasisKey(face,0);
522  LibUtilities::BasisKey tobasis1 = DetFaceBasisKey(face,1);
523 
524  // Number of quadrature points in face expansion.
525  int nq_face = tobasis0.GetNumPoints()*tobasis1.GetNumPoints();
526 
527  int vCoordDim = GetCoordim();
528  int i;
529 
530  m_faceNormals[face] = Array<OneD, Array<OneD, NekDouble> >(vCoordDim);
531  Array<OneD, Array<OneD, NekDouble> > &normal = m_faceNormals[face];
532  for (i = 0; i < vCoordDim; ++i)
533  {
534  normal[i] = Array<OneD, NekDouble>(nq_face);
535  }
536 
537  // Regular geometry case
538  if (type == SpatialDomains::eRegular ||
540  {
541  NekDouble fac;
542  // Set up normals
543  switch(face)
544  {
545  case 0:
546  {
547  for(i = 0; i < vCoordDim; ++i)
548  {
549  normal[i][0] = -df[3*i+2][0];
550  }
551  break;
552  }
553  case 1:
554  {
555  for(i = 0; i < vCoordDim; ++i)
556  {
557  normal[i][0] = -df[3*i+1][0];
558  }
559  break;
560  }
561  case 2:
562  {
563  for(i = 0; i < vCoordDim; ++i)
564  {
565  normal[i][0] = df[3*i][0]+df[3*i+2][0];
566  }
567  break;
568  }
569  case 3:
570  {
571  for(i = 0; i < vCoordDim; ++i)
572  {
573  normal[i][0] = df[3*i+1][0]+df[3*i+2][0];
574  }
575  break;
576  }
577  case 4:
578  {
579  for(i = 0; i < vCoordDim; ++i)
580  {
581  normal[i][0] = -df[3*i][0];
582  }
583  break;
584  }
585  default:
586  ASSERTL0(false,"face is out of range (face < 4)");
587  }
588 
589  // Normalise resulting vector.
590  fac = 0.0;
591  for(i = 0; i < vCoordDim; ++i)
592  {
593  fac += normal[i][0]*normal[i][0];
594  }
595  fac = 1.0/sqrt(fac);
596  for (i = 0; i < vCoordDim; ++i)
597  {
598  Vmath::Fill(nq_face,fac*normal[i][0],normal[i],1);
599  }
600  }
601  else
602  {
603  // Set up deformed normals.
604  int j, k;
605 
606  int nq0 = ptsKeys[0].GetNumPoints();
607  int nq1 = ptsKeys[1].GetNumPoints();
608  int nq2 = ptsKeys[2].GetNumPoints();
609  int nq01 = nq0*nq1;
610  int nqtot;
611 
612  // Determine number of quadrature points on the face.
613  if (face == 0)
614  {
615  nqtot = nq0*nq1;
616  }
617  else if (face == 1 || face == 3)
618  {
619  nqtot = nq0*nq2;
620  }
621  else
622  {
623  nqtot = nq1*nq2;
624  }
625 
626  LibUtilities::PointsKey points0;
627  LibUtilities::PointsKey points1;
628 
629  Array<OneD, NekDouble> faceJac(nqtot);
630  Array<OneD, NekDouble> normals(vCoordDim*nqtot,0.0);
631 
632  // Extract Jacobian along face and recover local derivatives
633  // (dx/dr) for polynomial interpolation by multiplying m_gmat by
634  // jacobian
635  switch(face)
636  {
637  case 0:
638  {
639  for(j = 0; j < nq01; ++j)
640  {
641  normals[j] = -df[2][j]*jac[j];
642  normals[nqtot+j] = -df[5][j]*jac[j];
643  normals[2*nqtot+j] = -df[8][j]*jac[j];
644  faceJac[j] = jac[j];
645  }
646 
647  points0 = ptsKeys[0];
648  points1 = ptsKeys[1];
649  break;
650  }
651 
652  case 1:
653  {
654  for (j = 0; j < nq0; ++j)
655  {
656  for(k = 0; k < nq2; ++k)
657  {
658  int tmp = j+nq01*k;
659  normals[j+k*nq0] =
660  -df[1][tmp]*jac[tmp];
661  normals[nqtot+j+k*nq0] =
662  -df[4][tmp]*jac[tmp];
663  normals[2*nqtot+j+k*nq0] =
664  -df[7][tmp]*jac[tmp];
665  faceJac[j+k*nq0] = jac[tmp];
666  }
667  }
668 
669  points0 = ptsKeys[0];
670  points1 = ptsKeys[2];
671  break;
672  }
673 
674  case 2:
675  {
676  for (j = 0; j < nq1; ++j)
677  {
678  for(k = 0; k < nq2; ++k)
679  {
680  int tmp = nq0-1+nq0*j+nq01*k;
681  normals[j+k*nq1] =
682  (df[0][tmp]+df[2][tmp])*jac[tmp];
683  normals[nqtot+j+k*nq1] =
684  (df[3][tmp]+df[5][tmp])*jac[tmp];
685  normals[2*nqtot+j+k*nq1] =
686  (df[6][tmp]+df[8][tmp])*jac[tmp];
687  faceJac[j+k*nq1] = jac[tmp];
688  }
689  }
690 
691  points0 = ptsKeys[1];
692  points1 = ptsKeys[2];
693  break;
694  }
695 
696  case 3:
697  {
698  for (j = 0; j < nq0; ++j)
699  {
700  for(k = 0; k < nq2; ++k)
701  {
702  int tmp = nq0*(nq1-1) + j + nq01*k;
703  normals[j+k*nq0] =
704  (df[1][tmp]+df[2][tmp])*jac[tmp];
705  normals[nqtot+j+k*nq0] =
706  (df[4][tmp]+df[5][tmp])*jac[tmp];
707  normals[2*nqtot+j+k*nq0] =
708  (df[7][tmp]+df[8][tmp])*jac[tmp];
709  faceJac[j+k*nq0] = jac[tmp];
710  }
711  }
712 
713  points0 = ptsKeys[0];
714  points1 = ptsKeys[2];
715  break;
716  }
717 
718  case 4:
719  {
720  for (j = 0; j < nq1; ++j)
721  {
722  for(k = 0; k < nq2; ++k)
723  {
724  int tmp = j*nq0+nq01*k;
725  normals[j+k*nq1] =
726  -df[0][tmp]*jac[tmp];
727  normals[nqtot+j+k*nq1] =
728  -df[3][tmp]*jac[tmp];
729  normals[2*nqtot+j+k*nq1] =
730  -df[6][tmp]*jac[tmp];
731  faceJac[j+k*nq1] = jac[tmp];
732  }
733  }
734 
735  points0 = ptsKeys[1];
736  points1 = ptsKeys[2];
737  break;
738  }
739 
740  default:
741  ASSERTL0(false,"face is out of range (face < 4)");
742  }
743 
744  Array<OneD, NekDouble> work (nq_face, 0.0);
745  // Interpolate Jacobian and invert
746  LibUtilities::Interp2D(points0, points1, faceJac,
747  tobasis0.GetPointsKey(),
748  tobasis1.GetPointsKey(),
749  work);
750  Vmath::Sdiv(nq_face, 1.0, &work[0], 1, &work[0], 1);
751 
752  // Interpolate normal and multiply by inverse Jacobian.
753  for(i = 0; i < vCoordDim; ++i)
754  {
755  LibUtilities::Interp2D(points0, points1,
756  &normals[i*nqtot],
757  tobasis0.GetPointsKey(),
758  tobasis1.GetPointsKey(),
759  &normal[i][0]);
760  Vmath::Vmul(nq_face,work,1,normal[i],1,normal[i],1);
761  }
762 
763  // Normalise to obtain unit normals.
764  Vmath::Zero(nq_face,work,1);
765  for(i = 0; i < GetCoordim(); ++i)
766  {
767  Vmath::Vvtvp(nq_face,normal[i],1,normal[i],1,work,1,work,1);
768  }
769 
770  Vmath::Vsqrt(nq_face,work,1,work,1);
771  Vmath::Sdiv (nq_face,1.0,work,1,work,1);
772 
773  for(i = 0; i < GetCoordim(); ++i)
774  {
775  Vmath::Vmul(nq_face,normal[i],1,work,1,normal[i],1);
776  }
777  }
778  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:220
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
Definition: Vmath.cpp:394
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:46
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:428
void Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/y.
Definition: Vmath.cpp:257
void Interp2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
this function interpolates a 2D function evaluated at the quadrature points of the 2D basis...
Definition: Interp.cpp:116
double NekDouble
std::map< int, NormalVector > m_faceNormals
SpatialDomains::GeometrySharedPtr GetGeom() const
Definition: Expansion.cpp:148
boost::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition: GeomFactors.h:62
Geometry is straight-sided with constant geometric factors.
const LibUtilities::BasisKey DetFaceBasisKey(const int i, const int k) const
Definition: StdExpansion.h:324
GeomType
Indicates the type of element geometry.
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:359
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:169
void Nektar::LocalRegions::PyrExp::v_ComputeLaplacianMetric ( )
protectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1070 of file PyrExp.cpp.

References Nektar::LocalRegions::Expansion::ComputeQuadratureMetric(), Nektar::SpatialDomains::eDeformed, Nektar::LocalRegions::eMetricLaplacian00, Nektar::LocalRegions::eMetricLaplacian01, Nektar::LocalRegions::eMetricLaplacian02, Nektar::LocalRegions::eMetricLaplacian11, Nektar::LocalRegions::eMetricLaplacian12, Nektar::LocalRegions::eMetricLaplacian22, Nektar::LocalRegions::eMetricQuadrature, Vmath::Fill(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::GetTotPoints(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_metricinfo, Nektar::LocalRegions::Expansion::m_metrics, Nektar::StdRegions::StdExpansion::MultiplyByQuadratureMetric(), Vmath::Svtsvtp(), Vmath::Svtvp(), Vmath::Vvtvp(), and Vmath::Vvtvvtp().

1071  {
1072  if (m_metrics.count(eMetricQuadrature) == 0)
1073  {
1075  }
1076 
1077  int i, j;
1078  const unsigned int nqtot = GetTotPoints();
1079  const unsigned int dim = 3;
1080  const MetricType m[3][3] = {
1084  };
1085 
1086  for (unsigned int i = 0; i < dim; ++i)
1087  {
1088  for (unsigned int j = i; j < dim; ++j)
1089  {
1090  m_metrics[m[i][j]] = Array<OneD, NekDouble>(nqtot);
1091  }
1092  }
1093 
1094  // Define shorthand synonyms for m_metrics storage
1095  Array<OneD,NekDouble> g0 (m_metrics[m[0][0]]);
1096  Array<OneD,NekDouble> g1 (m_metrics[m[1][1]]);
1097  Array<OneD,NekDouble> g2 (m_metrics[m[2][2]]);
1098  Array<OneD,NekDouble> g3 (m_metrics[m[0][1]]);
1099  Array<OneD,NekDouble> g4 (m_metrics[m[0][2]]);
1100  Array<OneD,NekDouble> g5 (m_metrics[m[1][2]]);
1101 
1102  // Allocate temporary storage
1103  Array<OneD,NekDouble> alloc(9*nqtot,0.0);
1104  Array<OneD,NekDouble> h0 (nqtot, alloc);
1105  Array<OneD,NekDouble> h1 (nqtot, alloc+ 1*nqtot);
1106  Array<OneD,NekDouble> h2 (nqtot, alloc+ 2*nqtot);
1107  Array<OneD,NekDouble> wsp1 (nqtot, alloc+ 3*nqtot);
1108  Array<OneD,NekDouble> wsp2 (nqtot, alloc+ 4*nqtot);
1109  Array<OneD,NekDouble> wsp3 (nqtot, alloc+ 5*nqtot);
1110  Array<OneD,NekDouble> wsp4 (nqtot, alloc+ 6*nqtot);
1111  Array<OneD,NekDouble> wsp5 (nqtot, alloc+ 7*nqtot);
1112  Array<OneD,NekDouble> wsp6 (nqtot, alloc+ 8*nqtot);
1113 
1114  const Array<TwoD, const NekDouble>& df =
1115  m_metricinfo->GetDerivFactors(GetPointsKeys());
1116  const Array<OneD, const NekDouble>& z0 = m_base[0]->GetZ();
1117  const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
1118  const Array<OneD, const NekDouble>& z2 = m_base[2]->GetZ();
1119  const unsigned int nquad0 = m_base[0]->GetNumPoints();
1120  const unsigned int nquad1 = m_base[1]->GetNumPoints();
1121  const unsigned int nquad2 = m_base[2]->GetNumPoints();
1122 
1123  // Populate collapsed coordinate arrays h0, h1 and h2.
1124  for(j = 0; j < nquad2; ++j)
1125  {
1126  for(i = 0; i < nquad1; ++i)
1127  {
1128  Vmath::Fill(nquad0, 2.0/(1.0-z2[j]), &h0[0]+i*nquad0 + j*nquad0*nquad1,1);
1129  Vmath::Fill(nquad0, 1.0/(1.0-z2[j]), &h1[0]+i*nquad0 + j*nquad0*nquad1,1);
1130  Vmath::Fill(nquad0, (1.0+z1[i])/(1.0-z2[j]), &h2[0]+i*nquad0 + j*nquad0*nquad1,1);
1131  }
1132  }
1133  for(i = 0; i < nquad0; i++)
1134  {
1135  Blas::Dscal(nquad1*nquad2, 1+z0[i], &h1[0]+i, nquad0);
1136  }
1137 
1138  // Step 3. Construct combined metric terms for physical space to
1139  // collapsed coordinate system.
1140  // Order of construction optimised to minimise temporary storage
1141  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1142  {
1143  // f_{1k}
1144  Vmath::Vvtvvtp(nqtot, &df[0][0], 1, &h0[0], 1, &df[2][0], 1, &h1[0], 1, &wsp1[0], 1);
1145  Vmath::Vvtvvtp(nqtot, &df[3][0], 1, &h0[0], 1, &df[5][0], 1, &h1[0], 1, &wsp2[0], 1);
1146  Vmath::Vvtvvtp(nqtot, &df[6][0], 1, &h0[0], 1, &df[8][0], 1, &h1[0], 1, &wsp3[0], 1);
1147 
1148  // g0
1149  Vmath::Vvtvvtp(nqtot, &wsp1[0], 1, &wsp1[0], 1, &wsp2[0], 1, &wsp2[0], 1, &g0[0], 1);
1150  Vmath::Vvtvp (nqtot, &wsp3[0], 1, &wsp3[0], 1, &g0[0], 1, &g0[0], 1);
1151 
1152  // g4
1153  Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &wsp1[0], 1, &df[5][0], 1, &wsp2[0], 1, &g4[0], 1);
1154  Vmath::Vvtvp (nqtot, &df[8][0], 1, &wsp3[0], 1, &g4[0], 1, &g4[0], 1);
1155 
1156  // f_{2k}
1157  Vmath::Vvtvvtp(nqtot, &df[1][0], 1, &h0[0], 1, &df[2][0], 1, &h2[0], 1, &wsp4[0], 1);
1158  Vmath::Vvtvvtp(nqtot, &df[4][0], 1, &h0[0], 1, &df[5][0], 1, &h2[0], 1, &wsp5[0], 1);
1159  Vmath::Vvtvvtp(nqtot, &df[7][0], 1, &h0[0], 1, &df[8][0], 1, &h2[0], 1, &wsp6[0], 1);
1160 
1161  // g1
1162  Vmath::Vvtvvtp(nqtot, &wsp4[0], 1, &wsp4[0], 1, &wsp5[0], 1, &wsp5[0], 1, &g1[0], 1);
1163  Vmath::Vvtvp (nqtot, &wsp6[0], 1, &wsp6[0], 1, &g1[0], 1, &g1[0], 1);
1164 
1165  // g3
1166  Vmath::Vvtvvtp(nqtot, &wsp1[0], 1, &wsp4[0], 1, &wsp2[0], 1, &wsp5[0], 1, &g3[0], 1);
1167  Vmath::Vvtvp (nqtot, &wsp3[0], 1, &wsp6[0], 1, &g3[0], 1, &g3[0], 1);
1168 
1169  // g5
1170  Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &wsp4[0], 1, &df[5][0], 1, &wsp5[0], 1, &g5[0], 1);
1171  Vmath::Vvtvp (nqtot, &df[8][0], 1, &wsp6[0], 1, &g5[0], 1, &g5[0], 1);
1172 
1173  // g2
1174  Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &df[2][0], 1, &df[5][0], 1, &df[5][0], 1, &g2[0], 1);
1175  Vmath::Vvtvp (nqtot, &df[8][0], 1, &df[8][0], 1, &g2[0], 1, &g2[0], 1);
1176  }
1177  else
1178  {
1179  // f_{1k}
1180  Vmath::Svtsvtp(nqtot, df[0][0], &h0[0], 1, df[2][0], &h1[0], 1, &wsp1[0], 1);
1181  Vmath::Svtsvtp(nqtot, df[3][0], &h0[0], 1, df[5][0], &h1[0], 1, &wsp2[0], 1);
1182  Vmath::Svtsvtp(nqtot, df[6][0], &h0[0], 1, df[8][0], &h1[0], 1, &wsp3[0], 1);
1183 
1184  // g0
1185  Vmath::Vvtvvtp(nqtot, &wsp1[0], 1, &wsp1[0], 1, &wsp2[0], 1, &wsp2[0], 1, &g0[0], 1);
1186  Vmath::Vvtvp (nqtot, &wsp3[0], 1, &wsp3[0], 1, &g0[0], 1, &g0[0], 1);
1187 
1188  // g4
1189  Vmath::Svtsvtp(nqtot, df[2][0], &wsp1[0], 1, df[5][0], &wsp2[0], 1, &g4[0], 1);
1190  Vmath::Svtvp (nqtot, df[8][0], &wsp3[0], 1, &g4[0], 1, &g4[0], 1);
1191 
1192  // f_{2k}
1193  Vmath::Svtsvtp(nqtot, df[1][0], &h0[0], 1, df[2][0], &h2[0], 1, &wsp4[0], 1);
1194  Vmath::Svtsvtp(nqtot, df[4][0], &h0[0], 1, df[5][0], &h2[0], 1, &wsp5[0], 1);
1195  Vmath::Svtsvtp(nqtot, df[7][0], &h0[0], 1, df[8][0], &h2[0], 1, &wsp6[0], 1);
1196 
1197  // g1
1198  Vmath::Vvtvvtp(nqtot, &wsp4[0], 1, &wsp4[0], 1, &wsp5[0], 1, &wsp5[0], 1, &g1[0], 1);
1199  Vmath::Vvtvp (nqtot, &wsp6[0], 1, &wsp6[0], 1, &g1[0], 1, &g1[0], 1);
1200 
1201  // g3
1202  Vmath::Vvtvvtp(nqtot, &wsp1[0], 1, &wsp4[0], 1, &wsp2[0], 1, &wsp5[0], 1, &g3[0], 1);
1203  Vmath::Vvtvp (nqtot, &wsp3[0], 1, &wsp6[0], 1, &g3[0], 1, &g3[0], 1);
1204 
1205  // g5
1206  Vmath::Svtsvtp(nqtot, df[2][0], &wsp4[0], 1, df[5][0], &wsp5[0], 1, &g5[0], 1);
1207  Vmath::Svtvp (nqtot, df[8][0], &wsp6[0], 1, &g5[0], 1, &g5[0], 1);
1208 
1209  // g2
1210  Vmath::Fill(nqtot, df[2][0]*df[2][0] + df[5][0]*df[5][0] + df[8][0]*df[8][0], &g2[0], 1);
1211  }
1212 
1213  for (unsigned int i = 0; i < dim; ++i)
1214  {
1215  for (unsigned int j = i; j < dim; ++j)
1216  {
1218  m_metrics[m[i][j]]);
1219 
1220  }
1221  }
1222  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:902
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:46
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:471
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:428
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:126
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:141
void Vvtvvtp(int n, const T *v, int incv, const T *w, int incw, const T *x, int incx, const T *y, int incy, T *z, int incz)
vvtvvtp (vector times vector plus vector times vector):
Definition: Vmath.cpp:523
void Svtsvtp(int n, const T alpha, const T *x, int incx, const T beta, const T *y, int incy, T *z, int incz)
vvtvvtp (scalar times vector plus scalar times vector):
Definition: Vmath.cpp:577
Array< OneD, LibUtilities::BasisSharedPtr > m_base
Geometry is curved or has non-constant factors.
DNekMatSharedPtr Nektar::LocalRegions::PyrExp::v_CreateStdMatrix ( const StdRegions::StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdPyrExp.

Definition at line 805 of file PyrExp.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and Nektar::StdRegions::StdExpansion::m_base.

806  {
807  LibUtilities::BasisKey bkey0 = m_base[0]->GetBasisKey();
808  LibUtilities::BasisKey bkey1 = m_base[1]->GetBasisKey();
809  LibUtilities::BasisKey bkey2 = m_base[2]->GetBasisKey();
811  MemoryManager<StdPyrExp>::AllocateSharedPtr(bkey0, bkey1, bkey2);
812 
813  return tmp->GetStdMatrix(mkey);
814  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< StdPyrExp > StdPyrExpSharedPtr
Definition: StdPyrExp.h:258
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::LocalRegions::PyrExp::v_DropLocStaticCondMatrix ( const MatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 826 of file PyrExp.cpp.

References m_staticCondMatrixManager.

827  {
828  m_staticCondMatrixManager.DeleteObject(mkey);
829  }
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: PyrExp.h:154
void Nektar::LocalRegions::PyrExp::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 (this)->m_coeffs.

Inputs:

  • inarray: array of physical quadrature points to be transformed

Outputs:

  • (this)->_coeffs: updated array of expansion coefficients.

Reimplemented from Nektar::StdRegions::StdPyrExp.

Definition at line 219 of file PyrExp.cpp.

References Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::StdRegions::eInvMass, Nektar::eWrapper, Nektar::StdRegions::StdExpansion::GetNcoeffs(), Nektar::StdRegions::StdExpansion::m_base, m_matrixManager, Nektar::StdRegions::StdExpansion::m_ncoeffs, v_IProductWRTBase(), and Vmath::Vcopy().

221  {
222  if(m_base[0]->Collocation() &&
223  m_base[1]->Collocation() &&
224  m_base[2]->Collocation())
225  {
226  Vmath::Vcopy(GetNcoeffs(),&inarray[0],1,&outarray[0],1);
227  }
228  else
229  {
230  v_IProductWRTBase(inarray,outarray);
231 
232  // get Mass matrix inverse
233  MatrixKey masskey(StdRegions::eInvMass,
234  DetShapeType(),*this);
235  DNekScalMatSharedPtr matsys = m_matrixManager[masskey];
236 
237  // copy inarray in case inarray == outarray
238  DNekVec in (m_ncoeffs,outarray);
239  DNekVec out(m_ncoeffs,outarray,eWrapper);
240 
241  out = (*matsys)*in;
242  }
243  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:454
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...
Definition: PyrExp.cpp:276
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
NekVector< NekDouble > DNekVec
Definition: NekTypeDefs.hpp:49
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:131
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1038
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: PyrExp.h:153
DNekMatSharedPtr Nektar::LocalRegions::PyrExp::v_GenMatrix ( const StdRegions::StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdPyrExp.

Definition at line 784 of file PyrExp.cpp.

References Nektar::StdRegions::eHybridDGHelmBndLam, Nektar::StdRegions::eHybridDGHelmholtz, Nektar::StdRegions::eHybridDGLamToQ0, Nektar::StdRegions::eHybridDGLamToQ1, Nektar::StdRegions::eHybridDGLamToQ2, Nektar::StdRegions::eHybridDGLamToU, Nektar::StdRegions::StdMatrixKey::GetMatrixType(), and Nektar::LocalRegions::Expansion3D::v_GenMatrix().

785  {
786  DNekMatSharedPtr returnval;
787 
788  switch(mkey.GetMatrixType())
789  {
796  returnval = Expansion3D::v_GenMatrix(mkey);
797  break;
798  default:
799  returnval = StdPyrExp::v_GenMatrix(mkey);
800  }
801 
802  return returnval;
803  }
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
virtual DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey)
void Nektar::LocalRegions::PyrExp::v_GetCoord ( const Array< OneD, const NekDouble > &  Lcoords,
Array< OneD, NekDouble > &  coords 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 316 of file PyrExp.cpp.

References ASSERTL1, and Nektar::LocalRegions::Expansion::m_geom.

318  {
319  int i;
320 
321  ASSERTL1(Lcoords[0] <= -1.0 && Lcoords[0] >= 1.0 &&
322  Lcoords[1] <= -1.0 && Lcoords[1] >= 1.0 &&
323  Lcoords[2] <= -1.0 && Lcoords[2] >= 1.0,
324  "Local coordinates are not in region [-1,1]");
325 
326  // m_geom->FillGeom(); // TODO: implement FillGeom()
327 
328  for(i = 0; i < m_geom->GetCoordim(); ++i)
329  {
330  coords[i] = m_geom->GetCoord(i,Lcoords);
331  }
332  }
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:125
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
int Nektar::LocalRegions::PyrExp::v_GetCoordim ( void  )
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 360 of file PyrExp.cpp.

References Nektar::LocalRegions::Expansion::m_geom.

361  {
362  return m_geom->GetCoordim();
363  }
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:125
void Nektar::LocalRegions::PyrExp::v_GetCoords ( Array< OneD, NekDouble > &  coords_1,
Array< OneD, NekDouble > &  coords_2,
Array< OneD, NekDouble > &  coords_3 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdPyrExp.

Definition at line 334 of file PyrExp.cpp.

References Nektar::LocalRegions::Expansion::v_GetCoords().

338  {
339  Expansion::v_GetCoords(coords_1, coords_2, coords_3);
340  }
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3)
Definition: Expansion.cpp:211
void Nektar::LocalRegions::PyrExp::v_GetFacePhysVals ( const int  face,
const StdRegions::StdExpansionSharedPtr FaceExp,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
StdRegions::Orientation  orient 
)
protectedvirtual

Definition at line 365 of file PyrExp.cpp.

References ASSERTL0, Nektar::StdRegions::eDir1BwdDir1_Dir2BwdDir2, Nektar::StdRegions::eDir1BwdDir1_Dir2FwdDir2, Nektar::StdRegions::eDir1BwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1BwdDir2_Dir2FwdDir1, Nektar::StdRegions::eDir1FwdDir1_Dir2BwdDir2, Nektar::StdRegions::eDir1FwdDir1_Dir2FwdDir2, Nektar::StdRegions::eDir1FwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1FwdDir2_Dir2FwdDir1, Nektar::StdRegions::eNoOrientation, Nektar::StdRegions::StdExpansion::GetFaceNumPoints(), Nektar::StdRegions::StdExpansion::GetForient(), Nektar::LibUtilities::Interp2D(), Nektar::StdRegions::StdExpansion::m_base, and Vmath::Vcopy().

371  {
372  int nq0 = m_base[0]->GetNumPoints();
373  int nq1 = m_base[1]->GetNumPoints();
374  int nq2 = m_base[2]->GetNumPoints();
375 
376  Array<OneD,NekDouble> o_tmp(GetFaceNumPoints(face));
377 
378  if (orient == StdRegions::eNoOrientation)
379  {
380  orient = GetForient(face);
381  }
382 
383  switch(face)
384  {
385  case 0:
387  {
388  //Directions A and B positive
389  Vmath::Vcopy(nq0*nq1,&(inarray[0]),1,&(o_tmp[0]),1);
390  }
391  else if(orient == StdRegions::eDir1BwdDir1_Dir2FwdDir2)
392  {
393  //Direction A negative and B positive
394  for (int j=0; j<nq1; j++)
395  {
396  Vmath::Vcopy(nq0,&(inarray[0])+(nq0-1)+j*nq0,-1,&(o_tmp[0])+(j*nq0),1);
397  }
398  }
399  else if(orient == StdRegions::eDir1FwdDir1_Dir2BwdDir2)
400  {
401  //Direction A positive and B negative
402  for (int j=0; j<nq1; j++)
403  {
404  Vmath::Vcopy(nq0,&(inarray[0])+nq0*(nq1-1-j),1,&(o_tmp[0])+(j*nq0),1);
405  }
406  }
407  else if(orient == StdRegions::eDir1BwdDir1_Dir2BwdDir2)
408  {
409  //Direction A negative and B negative
410  for(int j=0; j<nq1; j++)
411  {
412  Vmath::Vcopy(nq0,&(inarray[0])+(nq0*nq1-1-j*nq0),-1,&(o_tmp[0])+(j*nq0),1);
413  }
414  }
415  else if(orient == StdRegions::eDir1FwdDir2_Dir2FwdDir1)
416  {
417  //Transposed, Direction A and B positive
418  for (int i=0; i<nq0; i++)
419  {
420  Vmath::Vcopy(nq1,&(inarray[0])+i,nq0,&(o_tmp[0])+(i*nq1),1);
421  }
422  }
423  else if(orient == StdRegions::eDir1FwdDir2_Dir2BwdDir1)
424  {
425  //Transposed, Direction A positive and B negative
426  for (int i=0; i<nq0; i++)
427  {
428  Vmath::Vcopy(nq1,&(inarray[0])+(nq0-1-i),nq0,&(o_tmp[0])+(i*nq1),1);
429  }
430  }
431  else if(orient == StdRegions::eDir1BwdDir2_Dir2FwdDir1)
432  {
433  //Transposed, Direction A negative and B positive
434  for (int i=0; i<nq0; i++)
435  {
436  Vmath::Vcopy(nq1,&(inarray[0])+i+nq0*(nq1-1),-nq0,&(o_tmp[0])+(i*nq1),1);
437  }
438  }
439  else if(orient == StdRegions::eDir1BwdDir2_Dir2BwdDir1)
440  {
441  //Transposed, Direction A and B negative
442  for (int i=0; i<nq0; i++)
443  {
444  Vmath::Vcopy(nq1,&(inarray[0])+(nq0*nq1-1-i),-nq0,&(o_tmp[0])+(i*nq1),1);
445  }
446  }
447  LibUtilities::Interp2D(m_base[0]->GetPointsKey(), m_base[1]->GetPointsKey(), o_tmp,
448  FaceExp->GetBasis(0)->GetPointsKey(),FaceExp->GetBasis(1)->GetPointsKey(),outarray);
449  break;
450 
451  case 1:
452  {
453  for (int k = 0; k < nq2; k++)
454  {
455  Vmath::Vcopy(nq0,inarray.get()+nq0*nq1*k,1,outarray.get()+k*nq0,1);
456  }
457  LibUtilities::Interp2D(m_base[0]->GetPointsKey(), m_base[2]->GetPointsKey(), outarray.get(),
458  FaceExp->GetBasis(0)->GetPointsKey(),FaceExp->GetBasis(1)->GetPointsKey(),o_tmp.get());
459  break;
460  }
461 
462  case 2:
463  {
464  Vmath::Vcopy(nq1*nq2,inarray.get()+(nq0-1),nq0,outarray.get(),1);
465  LibUtilities::Interp2D(m_base[1]->GetPointsKey(), m_base[2]->GetPointsKey(), outarray.get(),
466  FaceExp->GetBasis(0)->GetPointsKey(),FaceExp->GetBasis(1)->GetPointsKey(),o_tmp.get());
467  break;
468  }
469 
470  case 3:
471  {
472  for (int k = 0; k < nq2; k++)
473  {
474  Vmath::Vcopy(nq0,inarray.get()+nq0*(nq1-1)+nq0*nq1*k,1,outarray.get()+(k*nq0),1);
475  }
476  LibUtilities::Interp2D(m_base[0]->GetPointsKey(), m_base[2]->GetPointsKey(), outarray.get(),
477  FaceExp->GetBasis(0)->GetPointsKey(),FaceExp->GetBasis(1)->GetPointsKey(),o_tmp.get());
478  }
479 
480  case 4:
481  {
482  Vmath::Vcopy(nq1*nq2,inarray.get(),nq0,outarray.get(),1);
483  LibUtilities::Interp2D(m_base[1]->GetPointsKey(), m_base[2]->GetPointsKey(), outarray.get(),
484  FaceExp->GetBasis(0)->GetPointsKey(),FaceExp->GetBasis(1)->GetPointsKey(),o_tmp.get());
485  break;
486  }
487 
488  default:
489  ASSERTL0(false,"face value (> 4) is out of range");
490  break;
491  }
492 
493  if (face > 0)
494  {
495  int fnq1 = FaceExp->GetNumPoints(0);
496  int fnq2 = FaceExp->GetNumPoints(1);
497 
498  if ((int)orient == 7)
499  {
500  for (int j = 0; j < fnq2; ++j)
501  {
502  Vmath::Vcopy(fnq1, o_tmp.get()+((j+1)*fnq1-1), -1, outarray.get()+j*fnq1, 1);
503  }
504  }
505  else
506  {
507  Vmath::Vcopy(fnq1*fnq2, o_tmp.get(), 1, outarray.get(), 1);
508  }
509  }
510  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
int GetFaceNumPoints(const int i) const
This function returns the number of quadrature points belonging to the i-th face. ...
Definition: StdExpansion.h:339
StdRegions::Orientation GetForient(int face)
Definition: StdExpansion.h:731
void Interp2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
this function interpolates a 2D function evaluated at the quadrature points of the 2D basis...
Definition: Interp.cpp:116
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1038
DNekScalMatSharedPtr Nektar::LocalRegions::PyrExp::v_GetLocMatrix ( const MatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 816 of file PyrExp.cpp.

References m_matrixManager.

817  {
818  return m_matrixManager[mkey];
819  }
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: PyrExp.h:153
DNekScalBlkMatSharedPtr Nektar::LocalRegions::PyrExp::v_GetLocStaticCondMatrix ( const MatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 821 of file PyrExp.cpp.

References m_staticCondMatrixManager.

822  {
823  return m_staticCondMatrixManager[mkey];
824  }
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: PyrExp.h:154
StdRegions::StdExpansionSharedPtr Nektar::LocalRegions::PyrExp::v_GetStdExp ( void  ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 304 of file PyrExp.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and Nektar::StdRegions::StdExpansion::m_base.

305  {
307  ::AllocateSharedPtr(m_base[0]->GetBasisKey(),
308  m_base[1]->GetBasisKey(),
309  m_base[2]->GetBasisKey());
310  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
Array< OneD, LibUtilities::BasisSharedPtr > m_base
NekDouble Nektar::LocalRegions::PyrExp::v_Integral ( const Array< OneD, const NekDouble > &  inarray)
protectedvirtual

Integrate the physical point list inarray over pyramidic region and return the value.

Inputs:

  • inarray: definition of function to be returned at quadrature point of expansion.

Outputs:

  • returns $\int^1_{-1}\int^1_{-1}\int^1_{-1} u(\bar \eta_1, \eta_2, \eta_3) J[i,j,k] d \bar \eta_1 d \eta_2 d \eta_3$
    $= \sum_{i=0}^{Q_1 - 1} \sum_{j=0}^{Q_2 - 1} \sum_{k=0}^{Q_3 - 1} u(\bar \eta_{1i}^{0,0}, \eta_{2j}^{0,0},\eta_{3k}^{2,0})w_{i}^{0,0} w_{j}^{0,0} \hat w_{k}^{2,0} $
    where $inarray[i,j, k] = u(\bar \eta_{1i},\eta_{2j}, \eta_{3k}) $,
    $\hat w_{k}^{2,0} = \frac {w^{2,0}} {2} $
    and $ J[i,j,k] $ is the Jacobian evaluated at the quadrature point.

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 110 of file PyrExp.cpp.

References Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_metricinfo, Vmath::Smul(), and Vmath::Vmul().

111  {
112  int nquad0 = m_base[0]->GetNumPoints();
113  int nquad1 = m_base[1]->GetNumPoints();
114  int nquad2 = m_base[2]->GetNumPoints();
115  Array<OneD, const NekDouble> jac = m_metricinfo->GetJac(GetPointsKeys());
116  Array<OneD, NekDouble> tmp(nquad0*nquad1*nquad2);
117 
118  // multiply inarray with Jacobian
119  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
120  {
121  Vmath::Vmul(nquad0*nquad1*nquad2,&jac[0],1,(NekDouble*)&inarray[0],1, &tmp[0],1);
122  }
123  else
124  {
125  Vmath::Smul(nquad0*nquad1*nquad2,(NekDouble) jac[0], (NekDouble*)&inarray[0],1,&tmp[0],1);
126  }
127 
128  // call StdPyrExp version;
129  return StdPyrExp::v_Integral(tmp);
130  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:126
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
double NekDouble
Array< OneD, LibUtilities::BasisSharedPtr > m_base
Geometry is curved or has non-constant factors.
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:169
void Nektar::LocalRegions::PyrExp::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} (\eta_{2j}) \psi_{pqr}^{c} (\eta_{3k}) w_i w_j w_k u(\bar \eta_{1,i} \eta_{2,j} \eta_{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(\eta_{2,j}) \sum_{k=0}^{nq_2} \psi_{pqr}^c u(\bar \eta_{1i},\eta_{2j},\eta_{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 (\eta_2) \psi_{pqr}^c (\eta_3) $
which can be implemented as
$f_{pqr} (\xi_{3k}) = \sum_{k=0}^{nq_3} \psi_{pqr}^c u(\bar \eta_{1i},\eta_{2j},\eta_{3k}) J_{i,j,k} = {\bf B_3 U} $
$ g_{pq} (\xi_{3k}) = \sum_{j=0}^{nq_1} \psi_{q}^a (\xi_{2j}) f_{pqr} (\xi_{3k}) = {\bf B_2 F} $
$ (\phi_{pqr}, u)_{\delta} = \sum_{k=0}^{nq_0} \psi_{p}^a (\xi_{3k}) g_{pq} (\xi_{3k}) = {\bf B_1 G} $

Reimplemented from Nektar::StdRegions::StdPyrExp.

Definition at line 276 of file PyrExp.cpp.

References Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_metricinfo, Vmath::Smul(), and Vmath::Vmul().

Referenced by v_FwdTrans().

279  {
280  int nquad0 = m_base[0]->GetNumPoints();
281  int nquad1 = m_base[1]->GetNumPoints();
282  int nquad2 = m_base[2]->GetNumPoints();
283  Array<OneD, const NekDouble> jac = m_metricinfo->GetJac(GetPointsKeys());
284  Array<OneD, NekDouble> tmp(nquad0*nquad1*nquad2);
285 
286  // multiply inarray with Jacobian
287  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
288  {
289  Vmath::Vmul(nquad0*nquad1*nquad2,&jac[0],1,(NekDouble*)&inarray[0],1,&tmp[0],1);
290  }
291  else
292  {
293  Vmath::Smul(nquad0*nquad1*nquad2,jac[0],(NekDouble*)&inarray[0],1,&tmp[0],1);
294  }
295 
296  StdPyrExp::v_IProductWRTBase(tmp,outarray);
297  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:126
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
double NekDouble
Array< OneD, LibUtilities::BasisSharedPtr > m_base
Geometry is curved or has non-constant factors.
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:169
void Nektar::LocalRegions::PyrExp::v_LaplacianMatrixOp_MatFree_Kernel ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
Array< OneD, NekDouble > &  wsp 
)
privatevirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1224 of file PyrExp.cpp.

References ASSERTL1, Nektar::LocalRegions::Expansion::ComputeLaplacianMetric(), Nektar::LocalRegions::eMetricLaplacian00, Nektar::LocalRegions::eMetricLaplacian01, Nektar::LocalRegions::eMetricLaplacian02, Nektar::LocalRegions::eMetricLaplacian11, Nektar::LocalRegions::eMetricLaplacian12, Nektar::LocalRegions::eMetricLaplacian22, Nektar::StdRegions::StdExpansion3D::IProductWRTBase_SumFacKernel(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_metrics, Nektar::StdRegions::StdExpansion::m_ncoeffs, Vmath::Vadd(), Vmath::Vvtvp(), and Vmath::Vvtvvtp().

1228  {
1229  // This implementation is only valid when there are no coefficients
1230  // associated to the Laplacian operator
1231  if (m_metrics.count(eMetricLaplacian00) == 0)
1232  {
1234  }
1235 
1236  int nquad0 = m_base[0]->GetNumPoints();
1237  int nquad1 = m_base[1]->GetNumPoints();
1238  int nq2 = m_base[2]->GetNumPoints();
1239  int nqtot = nquad0*nquad1*nq2;
1240 
1241  ASSERTL1(wsp.num_elements() >= 6*nqtot,
1242  "Insufficient workspace size.");
1243  ASSERTL1(m_ncoeffs <= nqtot,
1244  "Workspace not set up for ncoeffs > nqtot");
1245 
1246  const Array<OneD, const NekDouble>& base0 = m_base[0]->GetBdata();
1247  const Array<OneD, const NekDouble>& base1 = m_base[1]->GetBdata();
1248  const Array<OneD, const NekDouble>& base2 = m_base[2]->GetBdata();
1249  const Array<OneD, const NekDouble>& dbase0 = m_base[0]->GetDbdata();
1250  const Array<OneD, const NekDouble>& dbase1 = m_base[1]->GetDbdata();
1251  const Array<OneD, const NekDouble>& dbase2 = m_base[2]->GetDbdata();
1252  const Array<OneD, const NekDouble>& metric00 = m_metrics[eMetricLaplacian00];
1253  const Array<OneD, const NekDouble>& metric01 = m_metrics[eMetricLaplacian01];
1254  const Array<OneD, const NekDouble>& metric02 = m_metrics[eMetricLaplacian02];
1255  const Array<OneD, const NekDouble>& metric11 = m_metrics[eMetricLaplacian11];
1256  const Array<OneD, const NekDouble>& metric12 = m_metrics[eMetricLaplacian12];
1257  const Array<OneD, const NekDouble>& metric22 = m_metrics[eMetricLaplacian22];
1258 
1259  // Allocate temporary storage
1260  Array<OneD,NekDouble> wsp0 (2*nqtot, wsp);
1261  Array<OneD,NekDouble> wsp1 ( nqtot, wsp+1*nqtot);
1262  Array<OneD,NekDouble> wsp2 ( nqtot, wsp+2*nqtot);
1263  Array<OneD,NekDouble> wsp3 ( nqtot, wsp+3*nqtot);
1264  Array<OneD,NekDouble> wsp4 ( nqtot, wsp+4*nqtot);
1265  Array<OneD,NekDouble> wsp5 ( nqtot, wsp+5*nqtot);
1266 
1267  // LAPLACIAN MATRIX OPERATION
1268  // wsp1 = du_dxi1 = D_xi1 * inarray = D_xi1 * u
1269  // wsp2 = du_dxi2 = D_xi2 * inarray = D_xi2 * u
1270  // wsp2 = du_dxi3 = D_xi3 * inarray = D_xi3 * u
1271  StdExpansion3D::PhysTensorDeriv(inarray,wsp0,wsp1,wsp2);
1272 
1273  // wsp0 = k = g0 * wsp1 + g1 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
1274  // wsp2 = l = g1 * wsp1 + g2 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
1275  // where g0, g1 and g2 are the metric terms set up in the GeomFactors class
1276  // especially for this purpose
1277  Vmath::Vvtvvtp(nqtot,&metric00[0],1,&wsp0[0],1,&metric01[0],1,&wsp1[0],1,&wsp3[0],1);
1278  Vmath::Vvtvp (nqtot,&metric02[0],1,&wsp2[0],1,&wsp3[0],1,&wsp3[0],1);
1279  Vmath::Vvtvvtp(nqtot,&metric01[0],1,&wsp0[0],1,&metric11[0],1,&wsp1[0],1,&wsp4[0],1);
1280  Vmath::Vvtvp (nqtot,&metric12[0],1,&wsp2[0],1,&wsp4[0],1,&wsp4[0],1);
1281  Vmath::Vvtvvtp(nqtot,&metric02[0],1,&wsp0[0],1,&metric12[0],1,&wsp1[0],1,&wsp5[0],1);
1282  Vmath::Vvtvp (nqtot,&metric22[0],1,&wsp2[0],1,&wsp5[0],1,&wsp5[0],1);
1283 
1284  // outarray = m = (D_xi1 * B)^T * k
1285  // wsp1 = n = (D_xi2 * B)^T * l
1286  IProductWRTBase_SumFacKernel(dbase0,base1,base2,wsp3,outarray,wsp0,false,true,true);
1287  IProductWRTBase_SumFacKernel(base0,dbase1,base2,wsp4,wsp2, wsp0,true,false,true);
1288  Vmath::Vadd(m_ncoeffs,wsp2.get(),1,outarray.get(),1,outarray.get(),1);
1289  IProductWRTBase_SumFacKernel(base0,base1,dbase2,wsp5,wsp2, wsp0,true,true,false);
1290  Vmath::Vadd(m_ncoeffs,wsp2.get(),1,outarray.get(),1,outarray.get(),1);
1291  }
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:428
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 Vvtvvtp(int n, const T *v, int incv, const T *w, int incw, const T *x, int incx, const T *y, int incy, T *z, int incz)
vvtvvtp (vector times vector plus vector times vector):
Definition: Vmath.cpp:523
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:285
void Nektar::LocalRegions::PyrExp::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::StdPyrExp.

Definition at line 137 of file PyrExp.cpp.

References Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_metricinfo, Vmath::Smul(), Vmath::Vmul(), and Vmath::Vvtvp().

141  {
142  int nquad0 = m_base[0]->GetNumPoints();
143  int nquad1 = m_base[1]->GetNumPoints();
144  int nquad2 = m_base[2]->GetNumPoints();
145  Array<TwoD, const NekDouble> gmat =
146  m_metricinfo->GetDerivFactors(GetPointsKeys());
147  Array<OneD,NekDouble> diff0(nquad0*nquad1*nquad2);
148  Array<OneD,NekDouble> diff1(nquad0*nquad1*nquad2);
149  Array<OneD,NekDouble> diff2(nquad0*nquad1*nquad2);
150 
151  StdPyrExp::v_PhysDeriv(inarray, diff0, diff1, diff2);
152 
153  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
154  {
155  if(out_d0.num_elements())
156  {
157  Vmath::Vmul (nquad0*nquad1*nquad2,&gmat[0][0],1,&diff0[0],1, &out_d0[0], 1);
158  Vmath::Vvtvp (nquad0*nquad1*nquad2,&gmat[1][0],1,&diff1[0],1, &out_d0[0], 1,&out_d0[0],1);
159  Vmath::Vvtvp (nquad0*nquad1*nquad2,&gmat[2][0],1,&diff2[0],1, &out_d0[0], 1,&out_d0[0],1);
160  }
161 
162  if(out_d1.num_elements())
163  {
164  Vmath::Vmul (nquad0*nquad1*nquad2,&gmat[3][0],1,&diff0[0],1, &out_d1[0], 1);
165  Vmath::Vvtvp (nquad0*nquad1*nquad2,&gmat[4][0],1,&diff1[0],1, &out_d1[0], 1,&out_d1[0],1);
166  Vmath::Vvtvp (nquad0*nquad1*nquad2,&gmat[5][0],1,&diff2[0],1, &out_d1[0], 1,&out_d1[0],1);
167  }
168 
169  if(out_d2.num_elements())
170  {
171  Vmath::Vmul (nquad0*nquad1*nquad2,&gmat[6][0],1,&diff0[0],1, &out_d2[0], 1);
172  Vmath::Vvtvp (nquad0*nquad1*nquad2,&gmat[7][0],1,&diff1[0],1, &out_d2[0], 1, &out_d2[0],1);
173  Vmath::Vvtvp (nquad0*nquad1*nquad2,&gmat[8][0],1,&diff2[0],1, &out_d2[0], 1, &out_d2[0],1);
174  }
175  }
176  else // regular geometry
177  {
178  if(out_d0.num_elements())
179  {
180  Vmath::Smul (nquad0*nquad1*nquad2,gmat[0][0],&diff0[0],1, &out_d0[0], 1);
181  Blas::Daxpy (nquad0*nquad1*nquad2,gmat[1][0],&diff1[0],1, &out_d0[0], 1);
182  Blas::Daxpy (nquad0*nquad1*nquad2,gmat[2][0],&diff2[0],1, &out_d0[0], 1);
183  }
184 
185  if(out_d1.num_elements())
186  {
187  Vmath::Smul (nquad0*nquad1*nquad2,gmat[3][0],&diff0[0],1, &out_d1[0], 1);
188  Blas::Daxpy (nquad0*nquad1*nquad2,gmat[4][0],&diff1[0],1, &out_d1[0], 1);
189  Blas::Daxpy (nquad0*nquad1*nquad2,gmat[5][0],&diff2[0],1, &out_d1[0], 1);
190  }
191 
192  if(out_d2.num_elements())
193  {
194  Vmath::Smul (nquad0*nquad1*nquad2,gmat[6][0],&diff0[0],1, &out_d2[0], 1);
195  Blas::Daxpy (nquad0*nquad1*nquad2,gmat[7][0],&diff1[0],1, &out_d2[0], 1);
196  Blas::Daxpy (nquad0*nquad1*nquad2,gmat[8][0],&diff2[0],1, &out_d2[0], 1);
197  }
198  }
199  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:428
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:126
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
Array< OneD, LibUtilities::BasisSharedPtr > m_base
Geometry is curved or has non-constant factors.
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:169
NekDouble Nektar::LocalRegions::PyrExp::v_PhysEvaluate ( const Array< OneD, const NekDouble > &  coords,
const Array< OneD, const NekDouble > &  physvals 
)
protectedvirtual

This function evaluates the expansion at a single (arbitrary) point of the domain.

Based on the value of the expansion at the quadrature points, this function calculates the value of the expansion at an arbitrary single points (with coordinates $ \mathbf{x_c}$ given by the pointer coords). This operation, equivalent to

\[ u(\mathbf{x_c}) = \sum_p \phi_p(\mathbf{x_c}) \hat{u}_p \]

is evaluated using Lagrangian interpolants through the quadrature points:

\[ u(\mathbf{x_c}) = \sum_p h_p(\mathbf{x_c}) u_p\]

This function requires that the physical value array $\mathbf{u}$ (implemented as the attribute #phys) is set.

Parameters
coordsthe coordinates of the single point
Returns
returns the value of the expansion at the single point

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 342 of file PyrExp.cpp.

References ASSERTL0, and Nektar::LocalRegions::Expansion::m_geom.

344  {
345  Array<OneD,NekDouble> Lcoord(3);
346 
347  ASSERTL0(m_geom,"m_geom not defined");
348 
349  //TODO: check GetLocCoords()
350  m_geom->GetLocCoords(coord, Lcoord);
351 
352  return StdPyrExp::v_PhysEvaluate(Lcoord, physvals);
353  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:125

Member Data Documentation

LibUtilities::NekManager<MatrixKey, DNekScalMat, MatrixKey::opLess> Nektar::LocalRegions::PyrExp::m_matrixManager
private

Definition at line 153 of file PyrExp.h.

Referenced by CreateMatrix(), v_FwdTrans(), and v_GetLocMatrix().

LibUtilities::NekManager<MatrixKey, DNekScalBlkMat, MatrixKey::opLess> Nektar::LocalRegions::PyrExp::m_staticCondMatrixManager
private

Definition at line 154 of file PyrExp.h.

Referenced by v_DropLocStaticCondMatrix(), and v_GetLocStaticCondMatrix().