Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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::BasisSharedPtr
GetBasis (int dir) const
 This function gets the shared point to basis in the dir direction. More...
 
int GetNcoeffs (void) const
 This function returns the total number of coefficients used in the expansion. More...
 
int GetTotPoints () const
 This function returns the total number of quadrature points used in the element. More...
 
LibUtilities::BasisType GetBasisType (const int dir) const
 This function returns the type of basis used in the dir direction. More...
 
int GetBasisNumModes (const int dir) const
 This function returns the number of expansion modes in the dir direction. More...
 
int EvalBasisNumModesMax (void) const
 This function returns the maximum number of expansion modes over all local directions. More...
 
LibUtilities::PointsType GetPointsType (const int dir) const
 This function returns the type of quadrature points used in the dir direction. More...
 
int GetNumPoints (const int dir) const
 This function returns the number of quadrature points in the dir direction. More...
 
const Array< OneD, const
NekDouble > & 
GetPoints (const int dir) const
 This function returns a pointer to the array containing the quadrature points in dir direction. More...
 
int GetNverts () const
 This function returns the number of vertices of the expansion domain. More...
 
int GetNedges () const
 This function returns the number of edges of the expansion domain. More...
 
int GetEdgeNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th edge. More...
 
int GetTotalEdgeIntNcoeffs () const
 
int GetEdgeNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th edge. More...
 
int DetCartesianDirOfEdge (const int edge)
 
const LibUtilities::BasisKey DetEdgeBasisKey (const int i) const
 
const LibUtilities::BasisKey DetFaceBasisKey (const int i, const int k) const
 
int GetFaceNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th face. More...
 
int GetFaceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th face. More...
 
int GetFaceIntNcoeffs (const int i) const
 
int GetTotalFaceIntNcoeffs () const
 
int GetTraceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th edge/face. More...
 
LibUtilities::PointsKey GetFacePointsKey (const int i, const int j) const
 
int NumBndryCoeffs (void) const
 
int NumDGBndryCoeffs (void) const
 
LibUtilities::BasisType GetEdgeBasisType (const int i) const
 This function returns the type of expansion basis on the i-th edge. More...
 
const LibUtilities::PointsKey GetNodalPointsKey () const
 This function returns the type of expansion Nodal point type if defined. More...
 
int GetNfaces () const
 This function returns the number of faces of the expansion domain. More...
 
int GetNtrace () const
 Returns the number of trace elements connected to this element. More...
 
LibUtilities::ShapeType DetShapeType () const
 This function returns the shape of the expansion domain. More...
 
boost::shared_ptr< StdExpansionGetStdExp (void) const
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion ()
 
bool IsNodalNonTensorialExp ()
 
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space. More...
 
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space. More...
 
void FwdTrans_BndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble Integral (const Array< OneD, const NekDouble > &inarray)
 This function integrates the specified function over the domain. More...
 
void FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 This function fills the array outarray with the mode-th mode of the expansion. More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 this function calculates the inner product of a given function f with the different modes of the expansion More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check)
 
void IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
int GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void SetElmtId (const int id)
 Set the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2=NullNekDouble1DArray, Array< OneD, NekDouble > &coords_3=NullNekDouble1DArray)
 this function returns the physical coordinates of the quadrature points of the expansion More...
 
void GetCoord (const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
 given the coordinates of a point of the element in the local collapsed coordinate system, this function calculates the physical coordinates of the point More...
 
DNekMatSharedPtr GetStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr GetStdStaticCondMatrix (const StdMatrixKey &mkey)
 
IndexMapValuesSharedPtr GetIndexMap (const IndexMapKey &ikey)
 
const Array< OneD, const
NekDouble > & 
GetPhysNormals (void)
 
void SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
virtual void SetUpPhysNormals (const int edge)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
void DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
StdRegions::Orientation GetForient (int face)
 
StdRegions::Orientation GetEorient (int edge)
 
StdRegions::Orientation GetPorient (int point)
 
StdRegions::Orientation GetCartesianEorient (int edge)
 
void SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
 
void SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
int CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
void ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs)
 
NekDouble StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
int GetCoordim ()
 
void GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
void GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
int GetVertexMap (const int localVertexId, bool useCoeffPacking=false)
 
void GetEdgeInteriorMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
 
void GetFaceInteriorMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
 
void GetEdgeToElementMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int P=-1)
 
void GetFaceToElementMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int nummodesA=-1, int nummodesB=-1)
 
void GetEdgePhysVals (const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Extract the physical values along edge edge from inarray into outarray following the local edge orientation and point distribution defined by defined in EdgeExp. More...
 
void GetEdgePhysVals (const int edge, const boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetTracePhysVals (const int edge, const boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetVertexPhysVals (const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray)
 
void GetEdgeInterpVals (const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetEdgeQFactors (const int edge, Array< OneD, NekDouble > &outarray)
 Extract the metric factors to compute the contravariant fluxes along edge edge and stores them into outarray following the local edge orientation (i.e. anticlockwise convention). More...
 
void GetFacePhysVals (const int face, const boost::shared_ptr< StdExpansion > &FaceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=eNoOrientation)
 
void GetEdgePhysMap (const int edge, Array< OneD, int > &outarray)
 
void GetFacePhysMap (const int face, Array< OneD, int > &outarray)
 
void MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
DNekMatSharedPtr CreateGeneralMatrix (const StdMatrixKey &mkey)
 this function generates the mass matrix $\mathbf{M}[i][j] = \int \phi_i(\mathbf{x}) \phi_j(\mathbf{x}) d\mathbf{x}$ More...
 
void GeneralMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
DNekMatSharedPtr GenMatrix (const StdMatrixKey &mkey)
 
void PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void PhysDeriv_s (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds)
 
void PhysDeriv_n (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn)
 
void PhysDirectionalDeriv (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &outarray)
 
void StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void AddRobinMassMatrix (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
void AddRobinEdgeContribution (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, Array< OneD, NekDouble > &coeffs)
 
NekDouble PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta. More...
 
const boost::shared_ptr
< SpatialDomains::GeomFactors > & 
GetMetricInfo (void) const
 
virtual int v_GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
virtual const Array< OneD,
const NekDouble > & 
v_GetPhysNormals (void)
 
virtual void v_SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
virtual void v_SetUpPhysNormals (const int edge)
 
virtual void v_ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmode_offset, NekDouble *coeffs)
 Unpack data from input file assuming it comes from the same expansion type. More...
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
virtual StdRegions::Orientation v_GetEorient (int edge)
 
virtual StdRegions::Orientation v_GetCartesianEorient (int edge)
 
virtual StdRegions::Orientation v_GetPorient (int point)
 
NekDouble Linf (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete $ L_\infty$ error $ |\epsilon|_\infty = \max |u - u_{exact}|$ where $ u_{exact}$ is given by the array sol. More...
 
NekDouble L2 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete $ L_2$ error, $ | \epsilon |_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 dx \right]^{1/2} d\xi_1 $ where $ u_{exact}$ is given by the array sol. More...
 
NekDouble H1 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete $ H^1$ error, $ | \epsilon |^1_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 + \nabla(u - u_{exact})\cdot\nabla(u - u_{exact})\cdot dx \right]^{1/2} d\xi_1 $ where $ u_{exact}$ is given by the array sol. More...
 
const NormalVectorGetEdgeNormal (const int edge) const
 
void ComputeEdgeNormal (const int edge)
 
void NegateEdgeNormal (const int edge)
 
bool EdgeNormalNegated (const int edge)
 
void ComputeFaceNormal (const int face)
 
void NegateFaceNormal (const int face)
 
bool FaceNormalNegated (const int face)
 
void ComputeVertexNormal (const int vertex)
 
const NormalVectorGetFaceNormal (const int face) const
 
const NormalVectorGetVertexNormal (const int vertex) const
 
const NormalVectorGetSurfaceNormal (const int id) const
 
const LibUtilities::PointsKeyVector GetPointsKeys () const
 
Array< OneD, unsigned int > GetEdgeInverseBoundaryMap (int eid)
 
Array< OneD, unsigned int > GetFaceInverseBoundaryMap (int fid, StdRegions::Orientation faceOrient=eNoOrientation)
 
DNekMatSharedPtr BuildInverseTransformationMatrix (const DNekScalMatSharedPtr &m_transformationmatrix)
 
void PhysInterpToSimplexEquiSpaced (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int npset=-1)
 This function performs an interpolation from the physical space points provided at input into an array of equispaced points which are not the collapsed coordinate. So for a tetrahedron you will only get a tetrahedral number of values. More...
 
void GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true)
 This function provides the connectivity of local simplices (triangles or tets) to connect the equispaced data points provided by PhysInterpToSimplexEquiSpaced. More...
 
void EquiSpacedToCoeffs (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs a projection/interpolation from the equispaced points sometimes used in post-processing onto the coefficient space. More...
 
template<class T >
boost::shared_ptr< T > as ()
 
void IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
- 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
 
void ReOrientFacePhysMap (const int nvert, const StdRegions::Orientation orient, const int nq0, const int nq1, Array< OneD, int > &idmap)
 
void v_NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
- 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::GeomFactorsSharedPtr
v_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 ()
 
virtual void v_GetFacePhysMap (const int face, Array< OneD, int > &outarray)
 
void v_ComputeFaceNormal (const int face)
 
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)
 
virtual bool v_FaceNormalNegated (const int face)
 
virtual int v_GetTraceNcoeffs (const int i) const
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion
DNekMatSharedPtr CreateStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr CreateStdStaticCondMatrix (const StdMatrixKey &mkey)
 Create the static condensation of a matrix when using a boundary interior decomposition. More...
 
IndexMapValuesSharedPtr CreateIndexMap (const IndexMapKey &ikey)
 Create an IndexMap which contains mapping information linking any specific element shape with either its boundaries, edges, faces, verteces, etc. More...
 
void BwdTrans_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GeneralMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
 
void LaplacianMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp_MatFree (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void HelmholtzMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
 
virtual NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
- 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 void v_GetTracePhysVals (const int face, const StdRegions::StdExpansionSharedPtr &FaceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient)
 Returns the physical values at the quadrature points of a face Wrapper function to v_GetFacePhysVals. More...
 
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 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)
 
void ReOrientTriFacePhysMap (const StdRegions::Orientation orient, const int nq0, const int nq1, Array< OneD, int > &idmap)
 
void ReOrientQuadFacePhysMap (const StdRegions::Orientation orient, const int nq0, const int nq1, Array< OneD, int > &idmap)
 
- 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::opLess
m_matrixManager
 
LibUtilities::NekManager
< MatrixKey, DNekScalBlkMat,
MatrixKey::opLess
m_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::BasisSharedPtr
m_base
 
int m_elmt_id
 
int m_ncoeffs
 
LibUtilities::NekManager
< StdMatrixKey, DNekMat,
StdMatrixKey::opLess
m_stdMatrixManager
 
LibUtilities::NekManager
< StdMatrixKey, DNekBlkMat,
StdMatrixKey::opLess
m_stdStaticCondMatrixManager
 
LibUtilities::NekManager
< IndexMapKey, IndexMapValues,
IndexMapKey::opLess
m_IndexMapManager
 
- 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:910
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: PyrExp.h:152
DNekScalMatSharedPtr CreateMatrix(const MatrixKey &mkey)
Definition: PyrExp.cpp:784
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:151
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:152
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:151
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 784 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().

785  {
786  DNekScalMatSharedPtr returnval;
788 
789  ASSERTL2(m_metricinfo->GetGtype() != SpatialDomains::eNoGeomType,"Geometric information is not set up");
790 
791  switch(mkey.GetMatrixType())
792  {
793  case StdRegions::eMass:
794  {
795  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
796  {
797  NekDouble one = 1.0;
798  DNekMatSharedPtr mat = GenMatrix(mkey);
800  }
801  else
802  {
803  NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
804  DNekMatSharedPtr mat = GetStdMatrix(mkey);
806  }
807  }
808  break;
810  {
811  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
812  {
813  NekDouble one = 1.0;
814  StdRegions::StdMatrixKey masskey(StdRegions::eMass,DetShapeType(),
815  *this);
816  DNekMatSharedPtr mat = GenMatrix(masskey);
817  mat->Invert();
819  }
820  else
821  {
822  NekDouble fac = 1.0/(m_metricinfo->GetJac(ptsKeys))[0];
823  DNekMatSharedPtr mat = GetStdMatrix(mkey);
825  }
826  }
827  break;
829  {
830  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed ||
831  mkey.GetNVarCoeff() > 0 ||
832  mkey.ConstFactorExists(
834  {
835  NekDouble one = 1.0;
836  DNekMatSharedPtr mat = GenMatrix(mkey);
837 
839  }
840  else
841  {
842  MatrixKey lap00key(StdRegions::eLaplacian00,
843  mkey.GetShapeType(), *this);
844  MatrixKey lap01key(StdRegions::eLaplacian01,
845  mkey.GetShapeType(), *this);
846  MatrixKey lap02key(StdRegions::eLaplacian02,
847  mkey.GetShapeType(), *this);
848  MatrixKey lap11key(StdRegions::eLaplacian11,
849  mkey.GetShapeType(), *this);
850  MatrixKey lap12key(StdRegions::eLaplacian12,
851  mkey.GetShapeType(), *this);
852  MatrixKey lap22key(StdRegions::eLaplacian22,
853  mkey.GetShapeType(), *this);
854 
855  DNekMat &lap00 = *GetStdMatrix(lap00key);
856  DNekMat &lap01 = *GetStdMatrix(lap01key);
857  DNekMat &lap02 = *GetStdMatrix(lap02key);
858  DNekMat &lap11 = *GetStdMatrix(lap11key);
859  DNekMat &lap12 = *GetStdMatrix(lap12key);
860  DNekMat &lap22 = *GetStdMatrix(lap22key);
861 
862  NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
863  Array<TwoD, const NekDouble> gmat =
864  m_metricinfo->GetGmat(ptsKeys);
865 
866  int rows = lap00.GetRows();
867  int cols = lap00.GetColumns();
868 
870  ::AllocateSharedPtr(rows,cols);
871 
872  (*lap) = gmat[0][0]*lap00
873  + gmat[4][0]*lap11
874  + gmat[8][0]*lap22
875  + gmat[3][0]*(lap01 + Transpose(lap01))
876  + gmat[6][0]*(lap02 + Transpose(lap02))
877  + gmat[7][0]*(lap12 + Transpose(lap12));
878 
879  returnval = MemoryManager<DNekScalMat>
880  ::AllocateSharedPtr(jac, lap);
881  }
882  }
883  break;
885  {
886  NekDouble factor = mkey.GetConstFactor(StdRegions::eFactorLambda);
887  MatrixKey masskey(StdRegions::eMass, mkey.GetShapeType(), *this);
888  DNekScalMat &MassMat = *(this->m_matrixManager[masskey]);
889  MatrixKey lapkey(StdRegions::eLaplacian, mkey.GetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
890  DNekScalMat &LapMat = *(this->m_matrixManager[lapkey]);
891 
892  int rows = LapMat.GetRows();
893  int cols = LapMat.GetColumns();
894 
896 
897  (*helm) = LapMat + factor*MassMat;
898 
899  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(1.0, helm);
900  }
901  break;
902  default:
903  NEKERROR(ErrorUtil::efatal, "Matrix creation not defined");
904  break;
905  }
906 
907  return returnval;
908  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:470
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:700
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:151
DNekScalBlkMatSharedPtr Nektar::LocalRegions::PyrExp::CreateStaticCondMatrix ( const MatrixKey mkey)
protected

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

911  {
912  DNekScalBlkMatSharedPtr returnval;
913 
914  ASSERTL2(m_metricinfo->GetGtype() != SpatialDomains::eNoGeomType,"Geometric information is not set up");
915 
916  // set up block matrix system
917  unsigned int nbdry = NumBndryCoeffs();
918  unsigned int nint = (unsigned int)(m_ncoeffs - nbdry);
919  unsigned int exp_size[] = {nbdry, nint};
920  unsigned int nblks = 2;
921  returnval = MemoryManager<DNekScalBlkMat>::AllocateSharedPtr(nblks, nblks, exp_size, exp_size); //Really need a constructor which takes Arrays
922  NekDouble factor = 1.0;
923 
924  switch(mkey.GetMatrixType())
925  {
927  case StdRegions::eHelmholtz: // special case since Helmholtz not defined in StdRegions
928 
929  // use Deformed case for both regular and deformed geometries
930  factor = 1.0;
931  goto UseLocRegionsMatrix;
932  break;
933  default:
934  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
935  {
936  factor = 1.0;
937  goto UseLocRegionsMatrix;
938  }
939  else
940  {
942  factor = mat->Scale();
943  goto UseStdRegionsMatrix;
944  }
945  break;
946  UseStdRegionsMatrix:
947  {
948  NekDouble invfactor = 1.0/factor;
949  NekDouble one = 1.0;
952  DNekMatSharedPtr Asubmat;
953 
954  //TODO: check below
955  returnval->SetBlock(0,0,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(factor,Asubmat = mat->GetBlock(0,0)));
956  returnval->SetBlock(0,1,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,Asubmat = mat->GetBlock(0,1)));
957  returnval->SetBlock(1,0,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(factor,Asubmat = mat->GetBlock(1,0)));
958  returnval->SetBlock(1,1,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(invfactor,Asubmat = mat->GetBlock(1,1)));
959  }
960  break;
961  UseLocRegionsMatrix:
962  {
963  int i,j;
964  NekDouble invfactor = 1.0/factor;
965  NekDouble one = 1.0;
966  DNekScalMat &mat = *GetLocMatrix(mkey);
971 
972  Array<OneD,unsigned int> bmap(nbdry);
973  Array<OneD,unsigned int> imap(nint);
974  GetBoundaryMap(bmap);
975  GetInteriorMap(imap);
976 
977  for(i = 0; i < nbdry; ++i)
978  {
979  for(j = 0; j < nbdry; ++j)
980  {
981  (*A)(i,j) = mat(bmap[i],bmap[j]);
982  }
983 
984  for(j = 0; j < nint; ++j)
985  {
986  (*B)(i,j) = mat(bmap[i],imap[j]);
987  }
988  }
989 
990  for(i = 0; i < nint; ++i)
991  {
992  for(j = 0; j < nbdry; ++j)
993  {
994  (*C)(i,j) = mat(imap[i],bmap[j]);
995  }
996 
997  for(j = 0; j < nint; ++j)
998  {
999  (*D)(i,j) = mat(imap[i],imap[j]);
1000  }
1001  }
1002 
1003  // Calculate static condensed system
1004  if(nint)
1005  {
1006  D->Invert();
1007  (*B) = (*B)*(*D);
1008  (*A) = (*A) - (*B)*(*C);
1009  }
1010 
1011  DNekScalMatSharedPtr Atmp;
1012 
1013  returnval->SetBlock(0,0,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(factor,A));
1014  returnval->SetBlock(0,1,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,B));
1015  returnval->SetBlock(1,0,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(factor,C));
1016  returnval->SetBlock(1,1,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(invfactor,D));
1017 
1018  }
1019  }
1020  return returnval;
1021  }
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:705
boost::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:74
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:821
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:816
void Nektar::LocalRegions::PyrExp::v_ComputeFaceNormal ( const int  face)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

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

466  {
467  const SpatialDomains::GeomFactorsSharedPtr &geomFactors =
468  GetGeom()->GetMetricInfo();
470  SpatialDomains::GeomType type = geomFactors->GetGtype();
471  const Array<TwoD, const NekDouble> &df = geomFactors->GetDerivFactors(ptsKeys);
472  const Array<OneD, const NekDouble> &jac = geomFactors->GetJac(ptsKeys);
473 
474  LibUtilities::BasisKey tobasis0 = DetFaceBasisKey(face,0);
475  LibUtilities::BasisKey tobasis1 = DetFaceBasisKey(face,1);
476 
477  // Number of quadrature points in face expansion.
478  int nq_face = tobasis0.GetNumPoints()*tobasis1.GetNumPoints();
479 
480  int vCoordDim = GetCoordim();
481  int i;
482 
483  m_faceNormals[face] = Array<OneD, Array<OneD, NekDouble> >(vCoordDim);
484  Array<OneD, Array<OneD, NekDouble> > &normal = m_faceNormals[face];
485  for (i = 0; i < vCoordDim; ++i)
486  {
487  normal[i] = Array<OneD, NekDouble>(nq_face);
488  }
489 
490  // Regular geometry case
491  if (type == SpatialDomains::eRegular ||
493  {
494  NekDouble fac;
495  // Set up normals
496  switch(face)
497  {
498  case 0:
499  {
500  for(i = 0; i < vCoordDim; ++i)
501  {
502  normal[i][0] = -df[3*i+2][0];
503  }
504  break;
505  }
506  case 1:
507  {
508  for(i = 0; i < vCoordDim; ++i)
509  {
510  normal[i][0] = -df[3*i+1][0];
511  }
512  break;
513  }
514  case 2:
515  {
516  for(i = 0; i < vCoordDim; ++i)
517  {
518  normal[i][0] = df[3*i][0]+df[3*i+2][0];
519  }
520  break;
521  }
522  case 3:
523  {
524  for(i = 0; i < vCoordDim; ++i)
525  {
526  normal[i][0] = df[3*i+1][0]+df[3*i+2][0];
527  }
528  break;
529  }
530  case 4:
531  {
532  for(i = 0; i < vCoordDim; ++i)
533  {
534  normal[i][0] = -df[3*i][0];
535  }
536  break;
537  }
538  default:
539  ASSERTL0(false,"face is out of range (face < 4)");
540  }
541 
542  // Normalise resulting vector.
543  fac = 0.0;
544  for(i = 0; i < vCoordDim; ++i)
545  {
546  fac += normal[i][0]*normal[i][0];
547  }
548  fac = 1.0/sqrt(fac);
549  for (i = 0; i < vCoordDim; ++i)
550  {
551  Vmath::Fill(nq_face,fac*normal[i][0],normal[i],1);
552  }
553  }
554  else
555  {
556  // Set up deformed normals.
557  int j, k;
558 
559  int nq0 = ptsKeys[0].GetNumPoints();
560  int nq1 = ptsKeys[1].GetNumPoints();
561  int nq2 = ptsKeys[2].GetNumPoints();
562  int nq01 = nq0*nq1;
563  int nqtot;
564 
565  // Determine number of quadrature points on the face.
566  if (face == 0)
567  {
568  nqtot = nq0*nq1;
569  }
570  else if (face == 1 || face == 3)
571  {
572  nqtot = nq0*nq2;
573  }
574  else
575  {
576  nqtot = nq1*nq2;
577  }
578 
579  LibUtilities::PointsKey points0;
580  LibUtilities::PointsKey points1;
581 
582  Array<OneD, NekDouble> faceJac(nqtot);
583  Array<OneD, NekDouble> normals(vCoordDim*nqtot,0.0);
584 
585  // Extract Jacobian along face and recover local derivatives
586  // (dx/dr) for polynomial interpolation by multiplying m_gmat by
587  // jacobian
588  switch(face)
589  {
590  case 0:
591  {
592  for(j = 0; j < nq01; ++j)
593  {
594  normals[j] = -df[2][j]*jac[j];
595  normals[nqtot+j] = -df[5][j]*jac[j];
596  normals[2*nqtot+j] = -df[8][j]*jac[j];
597  faceJac[j] = jac[j];
598  }
599 
600  points0 = ptsKeys[0];
601  points1 = ptsKeys[1];
602  break;
603  }
604 
605  case 1:
606  {
607  for (j = 0; j < nq0; ++j)
608  {
609  for(k = 0; k < nq2; ++k)
610  {
611  int tmp = j+nq01*k;
612  normals[j+k*nq0] =
613  -df[1][tmp]*jac[tmp];
614  normals[nqtot+j+k*nq0] =
615  -df[4][tmp]*jac[tmp];
616  normals[2*nqtot+j+k*nq0] =
617  -df[7][tmp]*jac[tmp];
618  faceJac[j+k*nq0] = jac[tmp];
619  }
620  }
621 
622  points0 = ptsKeys[0];
623  points1 = ptsKeys[2];
624  break;
625  }
626 
627  case 2:
628  {
629  for (j = 0; j < nq1; ++j)
630  {
631  for(k = 0; k < nq2; ++k)
632  {
633  int tmp = nq0-1+nq0*j+nq01*k;
634  normals[j+k*nq1] =
635  (df[0][tmp]+df[2][tmp])*jac[tmp];
636  normals[nqtot+j+k*nq1] =
637  (df[3][tmp]+df[5][tmp])*jac[tmp];
638  normals[2*nqtot+j+k*nq1] =
639  (df[6][tmp]+df[8][tmp])*jac[tmp];
640  faceJac[j+k*nq1] = jac[tmp];
641  }
642  }
643 
644  points0 = ptsKeys[1];
645  points1 = ptsKeys[2];
646  break;
647  }
648 
649  case 3:
650  {
651  for (j = 0; j < nq0; ++j)
652  {
653  for(k = 0; k < nq2; ++k)
654  {
655  int tmp = nq0*(nq1-1) + j + nq01*k;
656  normals[j+k*nq0] =
657  (df[1][tmp]+df[2][tmp])*jac[tmp];
658  normals[nqtot+j+k*nq0] =
659  (df[4][tmp]+df[5][tmp])*jac[tmp];
660  normals[2*nqtot+j+k*nq0] =
661  (df[7][tmp]+df[8][tmp])*jac[tmp];
662  faceJac[j+k*nq0] = jac[tmp];
663  }
664  }
665 
666  points0 = ptsKeys[0];
667  points1 = ptsKeys[2];
668  break;
669  }
670 
671  case 4:
672  {
673  for (j = 0; j < nq1; ++j)
674  {
675  for(k = 0; k < nq2; ++k)
676  {
677  int tmp = j*nq0+nq01*k;
678  normals[j+k*nq1] =
679  -df[0][tmp]*jac[tmp];
680  normals[nqtot+j+k*nq1] =
681  -df[3][tmp]*jac[tmp];
682  normals[2*nqtot+j+k*nq1] =
683  -df[6][tmp]*jac[tmp];
684  faceJac[j+k*nq1] = jac[tmp];
685  }
686  }
687 
688  points0 = ptsKeys[1];
689  points1 = ptsKeys[2];
690  break;
691  }
692 
693  default:
694  ASSERTL0(false,"face is out of range (face < 4)");
695  }
696 
697  Array<OneD, NekDouble> work (nq_face, 0.0);
698  // Interpolate Jacobian and invert
699  LibUtilities::Interp2D(points0, points1, faceJac,
700  tobasis0.GetPointsKey(),
701  tobasis1.GetPointsKey(),
702  work);
703  Vmath::Sdiv(nq_face, 1.0, &work[0], 1, &work[0], 1);
704 
705  // Interpolate normal and multiply by inverse Jacobian.
706  for(i = 0; i < vCoordDim; ++i)
707  {
708  LibUtilities::Interp2D(points0, points1,
709  &normals[i*nqtot],
710  tobasis0.GetPointsKey(),
711  tobasis1.GetPointsKey(),
712  &normal[i][0]);
713  Vmath::Vmul(nq_face,work,1,normal[i],1,normal[i],1);
714  }
715 
716  // Normalise to obtain unit normals.
717  Vmath::Zero(nq_face,work,1);
718  for(i = 0; i < GetCoordim(); ++i)
719  {
720  Vmath::Vvtvp(nq_face,normal[i],1,normal[i],1,work,1,work,1);
721  }
722 
723  Vmath::Vsqrt(nq_face,work,1,work,1);
724  Vmath::Sdiv (nq_face,1.0,work,1,work,1);
725 
726  for(i = 0; i < GetCoordim(); ++i)
727  {
728  Vmath::Vmul(nq_face,normal[i],1,work,1,normal[i],1);
729  }
730  }
731  }
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 1023 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().

1024  {
1025  if (m_metrics.count(eMetricQuadrature) == 0)
1026  {
1028  }
1029 
1030  int i, j;
1031  const unsigned int nqtot = GetTotPoints();
1032  const unsigned int dim = 3;
1033  const MetricType m[3][3] = {
1037  };
1038 
1039  for (unsigned int i = 0; i < dim; ++i)
1040  {
1041  for (unsigned int j = i; j < dim; ++j)
1042  {
1043  m_metrics[m[i][j]] = Array<OneD, NekDouble>(nqtot);
1044  }
1045  }
1046 
1047  // Define shorthand synonyms for m_metrics storage
1048  Array<OneD,NekDouble> g0 (m_metrics[m[0][0]]);
1049  Array<OneD,NekDouble> g1 (m_metrics[m[1][1]]);
1050  Array<OneD,NekDouble> g2 (m_metrics[m[2][2]]);
1051  Array<OneD,NekDouble> g3 (m_metrics[m[0][1]]);
1052  Array<OneD,NekDouble> g4 (m_metrics[m[0][2]]);
1053  Array<OneD,NekDouble> g5 (m_metrics[m[1][2]]);
1054 
1055  // Allocate temporary storage
1056  Array<OneD,NekDouble> alloc(9*nqtot,0.0);
1057  Array<OneD,NekDouble> h0 (nqtot, alloc);
1058  Array<OneD,NekDouble> h1 (nqtot, alloc+ 1*nqtot);
1059  Array<OneD,NekDouble> h2 (nqtot, alloc+ 2*nqtot);
1060  Array<OneD,NekDouble> wsp1 (nqtot, alloc+ 3*nqtot);
1061  Array<OneD,NekDouble> wsp2 (nqtot, alloc+ 4*nqtot);
1062  Array<OneD,NekDouble> wsp3 (nqtot, alloc+ 5*nqtot);
1063  Array<OneD,NekDouble> wsp4 (nqtot, alloc+ 6*nqtot);
1064  Array<OneD,NekDouble> wsp5 (nqtot, alloc+ 7*nqtot);
1065  Array<OneD,NekDouble> wsp6 (nqtot, alloc+ 8*nqtot);
1066 
1067  const Array<TwoD, const NekDouble>& df =
1068  m_metricinfo->GetDerivFactors(GetPointsKeys());
1069  const Array<OneD, const NekDouble>& z0 = m_base[0]->GetZ();
1070  const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
1071  const Array<OneD, const NekDouble>& z2 = m_base[2]->GetZ();
1072  const unsigned int nquad0 = m_base[0]->GetNumPoints();
1073  const unsigned int nquad1 = m_base[1]->GetNumPoints();
1074  const unsigned int nquad2 = m_base[2]->GetNumPoints();
1075 
1076  // Populate collapsed coordinate arrays h0, h1 and h2.
1077  for(j = 0; j < nquad2; ++j)
1078  {
1079  for(i = 0; i < nquad1; ++i)
1080  {
1081  Vmath::Fill(nquad0, 2.0/(1.0-z2[j]), &h0[0]+i*nquad0 + j*nquad0*nquad1,1);
1082  Vmath::Fill(nquad0, 1.0/(1.0-z2[j]), &h1[0]+i*nquad0 + j*nquad0*nquad1,1);
1083  Vmath::Fill(nquad0, (1.0+z1[i])/(1.0-z2[j]), &h2[0]+i*nquad0 + j*nquad0*nquad1,1);
1084  }
1085  }
1086  for(i = 0; i < nquad0; i++)
1087  {
1088  Blas::Dscal(nquad1*nquad2, 1+z0[i], &h1[0]+i, nquad0);
1089  }
1090 
1091  // Step 3. Construct combined metric terms for physical space to
1092  // collapsed coordinate system.
1093  // Order of construction optimised to minimise temporary storage
1094  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1095  {
1096  // f_{1k}
1097  Vmath::Vvtvvtp(nqtot, &df[0][0], 1, &h0[0], 1, &df[2][0], 1, &h1[0], 1, &wsp1[0], 1);
1098  Vmath::Vvtvvtp(nqtot, &df[3][0], 1, &h0[0], 1, &df[5][0], 1, &h1[0], 1, &wsp2[0], 1);
1099  Vmath::Vvtvvtp(nqtot, &df[6][0], 1, &h0[0], 1, &df[8][0], 1, &h1[0], 1, &wsp3[0], 1);
1100 
1101  // g0
1102  Vmath::Vvtvvtp(nqtot, &wsp1[0], 1, &wsp1[0], 1, &wsp2[0], 1, &wsp2[0], 1, &g0[0], 1);
1103  Vmath::Vvtvp (nqtot, &wsp3[0], 1, &wsp3[0], 1, &g0[0], 1, &g0[0], 1);
1104 
1105  // g4
1106  Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &wsp1[0], 1, &df[5][0], 1, &wsp2[0], 1, &g4[0], 1);
1107  Vmath::Vvtvp (nqtot, &df[8][0], 1, &wsp3[0], 1, &g4[0], 1, &g4[0], 1);
1108 
1109  // f_{2k}
1110  Vmath::Vvtvvtp(nqtot, &df[1][0], 1, &h0[0], 1, &df[2][0], 1, &h2[0], 1, &wsp4[0], 1);
1111  Vmath::Vvtvvtp(nqtot, &df[4][0], 1, &h0[0], 1, &df[5][0], 1, &h2[0], 1, &wsp5[0], 1);
1112  Vmath::Vvtvvtp(nqtot, &df[7][0], 1, &h0[0], 1, &df[8][0], 1, &h2[0], 1, &wsp6[0], 1);
1113 
1114  // g1
1115  Vmath::Vvtvvtp(nqtot, &wsp4[0], 1, &wsp4[0], 1, &wsp5[0], 1, &wsp5[0], 1, &g1[0], 1);
1116  Vmath::Vvtvp (nqtot, &wsp6[0], 1, &wsp6[0], 1, &g1[0], 1, &g1[0], 1);
1117 
1118  // g3
1119  Vmath::Vvtvvtp(nqtot, &wsp1[0], 1, &wsp4[0], 1, &wsp2[0], 1, &wsp5[0], 1, &g3[0], 1);
1120  Vmath::Vvtvp (nqtot, &wsp3[0], 1, &wsp6[0], 1, &g3[0], 1, &g3[0], 1);
1121 
1122  // g5
1123  Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &wsp4[0], 1, &df[5][0], 1, &wsp5[0], 1, &g5[0], 1);
1124  Vmath::Vvtvp (nqtot, &df[8][0], 1, &wsp6[0], 1, &g5[0], 1, &g5[0], 1);
1125 
1126  // g2
1127  Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &df[2][0], 1, &df[5][0], 1, &df[5][0], 1, &g2[0], 1);
1128  Vmath::Vvtvp (nqtot, &df[8][0], 1, &df[8][0], 1, &g2[0], 1, &g2[0], 1);
1129  }
1130  else
1131  {
1132  // f_{1k}
1133  Vmath::Svtsvtp(nqtot, df[0][0], &h0[0], 1, df[2][0], &h1[0], 1, &wsp1[0], 1);
1134  Vmath::Svtsvtp(nqtot, df[3][0], &h0[0], 1, df[5][0], &h1[0], 1, &wsp2[0], 1);
1135  Vmath::Svtsvtp(nqtot, df[6][0], &h0[0], 1, df[8][0], &h1[0], 1, &wsp3[0], 1);
1136 
1137  // g0
1138  Vmath::Vvtvvtp(nqtot, &wsp1[0], 1, &wsp1[0], 1, &wsp2[0], 1, &wsp2[0], 1, &g0[0], 1);
1139  Vmath::Vvtvp (nqtot, &wsp3[0], 1, &wsp3[0], 1, &g0[0], 1, &g0[0], 1);
1140 
1141  // g4
1142  Vmath::Svtsvtp(nqtot, df[2][0], &wsp1[0], 1, df[5][0], &wsp2[0], 1, &g4[0], 1);
1143  Vmath::Svtvp (nqtot, df[8][0], &wsp3[0], 1, &g4[0], 1, &g4[0], 1);
1144 
1145  // f_{2k}
1146  Vmath::Svtsvtp(nqtot, df[1][0], &h0[0], 1, df[2][0], &h2[0], 1, &wsp4[0], 1);
1147  Vmath::Svtsvtp(nqtot, df[4][0], &h0[0], 1, df[5][0], &h2[0], 1, &wsp5[0], 1);
1148  Vmath::Svtsvtp(nqtot, df[7][0], &h0[0], 1, df[8][0], &h2[0], 1, &wsp6[0], 1);
1149 
1150  // g1
1151  Vmath::Vvtvvtp(nqtot, &wsp4[0], 1, &wsp4[0], 1, &wsp5[0], 1, &wsp5[0], 1, &g1[0], 1);
1152  Vmath::Vvtvp (nqtot, &wsp6[0], 1, &wsp6[0], 1, &g1[0], 1, &g1[0], 1);
1153 
1154  // g3
1155  Vmath::Vvtvvtp(nqtot, &wsp1[0], 1, &wsp4[0], 1, &wsp2[0], 1, &wsp5[0], 1, &g3[0], 1);
1156  Vmath::Vvtvp (nqtot, &wsp3[0], 1, &wsp6[0], 1, &g3[0], 1, &g3[0], 1);
1157 
1158  // g5
1159  Vmath::Svtsvtp(nqtot, df[2][0], &wsp4[0], 1, df[5][0], &wsp5[0], 1, &g5[0], 1);
1160  Vmath::Svtvp (nqtot, df[8][0], &wsp6[0], 1, &g5[0], 1, &g5[0], 1);
1161 
1162  // g2
1163  Vmath::Fill(nqtot, df[2][0]*df[2][0] + df[5][0]*df[5][0] + df[8][0]*df[8][0], &g2[0], 1);
1164  }
1165 
1166  for (unsigned int i = 0; i < dim; ++i)
1167  {
1168  for (unsigned int j = i; j < dim; ++j)
1169  {
1171  m_metrics[m[i][j]]);
1172 
1173  }
1174  }
1175  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:942
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 758 of file PyrExp.cpp.

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

759  {
760  LibUtilities::BasisKey bkey0 = m_base[0]->GetBasisKey();
761  LibUtilities::BasisKey bkey1 = m_base[1]->GetBasisKey();
762  LibUtilities::BasisKey bkey2 = m_base[2]->GetBasisKey();
764  MemoryManager<StdPyrExp>::AllocateSharedPtr(bkey0, bkey1, bkey2);
765 
766  return tmp->GetStdMatrix(mkey);
767  }
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 779 of file PyrExp.cpp.

References m_staticCondMatrixManager.

780  {
781  m_staticCondMatrixManager.DeleteObject(mkey);
782  }
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: PyrExp.h:152
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:470
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Calculate the inner product of inarray with respect to the basis B=base0*base1*base2 and put into out...
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:1047
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: PyrExp.h:151
DNekMatSharedPtr Nektar::LocalRegions::PyrExp::v_GenMatrix ( const StdRegions::StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdPyrExp.

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

738  {
739  DNekMatSharedPtr returnval;
740 
741  switch(mkey.GetMatrixType())
742  {
749  returnval = Expansion3D::v_GenMatrix(mkey);
750  break;
751  default:
752  returnval = StdPyrExp::v_GenMatrix(mkey);
753  }
754 
755  return returnval;
756  }
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_GetFacePhysMap ( const int  face,
Array< OneD, int > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 365 of file PyrExp.cpp.

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

367  {
368  int nquad0 = m_base[0]->GetNumPoints();
369  int nquad1 = m_base[1]->GetNumPoints();
370  int nquad2 = m_base[2]->GetNumPoints();
371 
372  int nq0 = 0;
373  int nq1 = 0;
374 
375  switch(face)
376  {
377  case 0:
378  nq0 = nquad0;
379  nq1 = nquad1;
380  if(outarray.num_elements()!=nq0*nq1)
381  {
382  outarray = Array<OneD, int>(nq0*nq1);
383  }
384 
385  //Directions A and B positive
386  for(int i = 0; i < nquad0*nquad1; ++i)
387  {
388  outarray[i] = i;
389  }
390 
391  break;
392  case 1:
393  nq0 = nquad0;
394  nq1 = nquad2;
395  if(outarray.num_elements()!=nq0*nq1)
396  {
397  outarray = Array<OneD, int>(nq0*nq1);
398  }
399 
400  //Direction A and B positive
401  for (int k=0; k<nquad2; k++)
402  {
403  for(int i = 0; i < nquad0; ++i)
404  {
405  outarray[k*nquad0+i] = (nquad0*nquad1*k)+i;
406  }
407  }
408 
409  break;
410  case 2:
411  nq0 = nquad1;
412  nq1 = nquad2;
413  if(outarray.num_elements()!=nq0*nq1)
414  {
415  outarray = Array<OneD, int>(nq0*nq1);
416  }
417 
418  //Directions A and B positive
419  for(int j = 0; j < nquad1*nquad2; ++j)
420  {
421  outarray[j] = nquad0-1 + j*nquad0;
422 
423  }
424  break;
425  case 3:
426 
427  nq0 = nquad0;
428  nq1 = nquad2;
429  if(outarray.num_elements()!=nq0*nq1)
430  {
431  outarray = Array<OneD, int>(nq0*nq1);
432  }
433 
434  //Direction A and B positive
435  for (int k=0; k<nquad2; k++)
436  {
437  for(int i = 0; i < nquad0; ++i)
438  {
439  outarray[k*nquad0+i] = nquad0*(nquad1-1) + (nquad0*nquad1*k)+i;
440  }
441  }
442 
443  case 4:
444  nq0 = nquad1;
445  nq1 = nquad2;
446 
447  if(outarray.num_elements()!=nq0*nq1)
448  {
449  outarray = Array<OneD, int>(nq0*nq1);
450  }
451 
452  //Directions A and B positive
453  for(int j = 0; j < nquad1*nquad2; ++j)
454  {
455  outarray[j] = j*nquad0;
456 
457  }
458  break;
459  default:
460  ASSERTL0(false,"face value (> 4) is out of range");
461  break;
462  }
463  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
Array< OneD, LibUtilities::BasisSharedPtr > m_base
DNekScalMatSharedPtr Nektar::LocalRegions::PyrExp::v_GetLocMatrix ( const MatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 769 of file PyrExp.cpp.

References m_matrixManager.

770  {
771  return m_matrixManager[mkey];
772  }
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: PyrExp.h:151
DNekScalBlkMatSharedPtr Nektar::LocalRegions::PyrExp::v_GetLocStaticCondMatrix ( const MatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 774 of file PyrExp.cpp.

References m_staticCondMatrixManager.

775  {
776  return m_staticCondMatrixManager[mkey];
777  }
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: PyrExp.h:152
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 1177 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().

1181  {
1182  // This implementation is only valid when there are no coefficients
1183  // associated to the Laplacian operator
1184  if (m_metrics.count(eMetricLaplacian00) == 0)
1185  {
1187  }
1188 
1189  int nquad0 = m_base[0]->GetNumPoints();
1190  int nquad1 = m_base[1]->GetNumPoints();
1191  int nq2 = m_base[2]->GetNumPoints();
1192  int nqtot = nquad0*nquad1*nq2;
1193 
1194  ASSERTL1(wsp.num_elements() >= 6*nqtot,
1195  "Insufficient workspace size.");
1196  ASSERTL1(m_ncoeffs <= nqtot,
1197  "Workspace not set up for ncoeffs > nqtot");
1198 
1199  const Array<OneD, const NekDouble>& base0 = m_base[0]->GetBdata();
1200  const Array<OneD, const NekDouble>& base1 = m_base[1]->GetBdata();
1201  const Array<OneD, const NekDouble>& base2 = m_base[2]->GetBdata();
1202  const Array<OneD, const NekDouble>& dbase0 = m_base[0]->GetDbdata();
1203  const Array<OneD, const NekDouble>& dbase1 = m_base[1]->GetDbdata();
1204  const Array<OneD, const NekDouble>& dbase2 = m_base[2]->GetDbdata();
1205  const Array<OneD, const NekDouble>& metric00 = m_metrics[eMetricLaplacian00];
1206  const Array<OneD, const NekDouble>& metric01 = m_metrics[eMetricLaplacian01];
1207  const Array<OneD, const NekDouble>& metric02 = m_metrics[eMetricLaplacian02];
1208  const Array<OneD, const NekDouble>& metric11 = m_metrics[eMetricLaplacian11];
1209  const Array<OneD, const NekDouble>& metric12 = m_metrics[eMetricLaplacian12];
1210  const Array<OneD, const NekDouble>& metric22 = m_metrics[eMetricLaplacian22];
1211 
1212  // Allocate temporary storage
1213  Array<OneD,NekDouble> wsp0 (2*nqtot, wsp);
1214  Array<OneD,NekDouble> wsp1 ( nqtot, wsp+1*nqtot);
1215  Array<OneD,NekDouble> wsp2 ( nqtot, wsp+2*nqtot);
1216  Array<OneD,NekDouble> wsp3 ( nqtot, wsp+3*nqtot);
1217  Array<OneD,NekDouble> wsp4 ( nqtot, wsp+4*nqtot);
1218  Array<OneD,NekDouble> wsp5 ( nqtot, wsp+5*nqtot);
1219 
1220  // LAPLACIAN MATRIX OPERATION
1221  // wsp1 = du_dxi1 = D_xi1 * inarray = D_xi1 * u
1222  // wsp2 = du_dxi2 = D_xi2 * inarray = D_xi2 * u
1223  // wsp2 = du_dxi3 = D_xi3 * inarray = D_xi3 * u
1224  StdExpansion3D::PhysTensorDeriv(inarray,wsp0,wsp1,wsp2);
1225 
1226  // wsp0 = k = g0 * wsp1 + g1 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
1227  // wsp2 = l = g1 * wsp1 + g2 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
1228  // where g0, g1 and g2 are the metric terms set up in the GeomFactors class
1229  // especially for this purpose
1230  Vmath::Vvtvvtp(nqtot,&metric00[0],1,&wsp0[0],1,&metric01[0],1,&wsp1[0],1,&wsp3[0],1);
1231  Vmath::Vvtvp (nqtot,&metric02[0],1,&wsp2[0],1,&wsp3[0],1,&wsp3[0],1);
1232  Vmath::Vvtvvtp(nqtot,&metric01[0],1,&wsp0[0],1,&metric11[0],1,&wsp1[0],1,&wsp4[0],1);
1233  Vmath::Vvtvp (nqtot,&metric12[0],1,&wsp2[0],1,&wsp4[0],1,&wsp4[0],1);
1234  Vmath::Vvtvvtp(nqtot,&metric02[0],1,&wsp0[0],1,&metric12[0],1,&wsp1[0],1,&wsp5[0],1);
1235  Vmath::Vvtvp (nqtot,&metric22[0],1,&wsp2[0],1,&wsp5[0],1,&wsp5[0],1);
1236 
1237  // outarray = m = (D_xi1 * B)^T * k
1238  // wsp1 = n = (D_xi2 * B)^T * l
1239  IProductWRTBase_SumFacKernel(dbase0,base1,base2,wsp3,outarray,wsp0,false,true,true);
1240  IProductWRTBase_SumFacKernel(base0,dbase1,base2,wsp4,wsp2, wsp0,true,false,true);
1241  Vmath::Vadd(m_ncoeffs,wsp2.get(),1,outarray.get(),1,outarray.get(),1);
1242  IProductWRTBase_SumFacKernel(base0,base1,dbase2,wsp5,wsp2, wsp0,true,true,false);
1243  Vmath::Vadd(m_ncoeffs,wsp2.get(),1,outarray.get(),1,outarray.get(),1);
1244  }
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 151 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 152 of file PyrExp.h.

Referenced by v_DropLocStaticCondMatrix(), and v_GetLocStaticCondMatrix().