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

#include <PyrExp.h>

Inheritance diagram for Nektar::LocalRegions::PyrExp:
[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 ()
 
- 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)
 
int GetNedges () const
 return the number of edges in 3D expansion More...
 
int GetEdgeNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th edge. More...
 
void GetEdgeInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards)
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion
 StdExpansion ()
 Default Constructor. More...
 
 StdExpansion (const int numcoeffs, const int numbases, const LibUtilities::BasisKey &Ba=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bb=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bc=LibUtilities::NullBasisKey)
 Constructor. More...
 
 StdExpansion (const StdExpansion &T)
 Copy Constructor. More...
 
virtual ~StdExpansion ()
 Destructor. More...
 
int GetNumBases () const
 This function returns the number of 1D bases used in the expansion. More...
 
const Array< OneD, const LibUtilities::BasisSharedPtr > & GetBase () const
 This function gets the shared point to basis. More...
 
const LibUtilities::BasisSharedPtrGetBasis (int dir) const
 This function gets the shared point to basis in the dir direction. More...
 
int GetNcoeffs (void) const
 This function returns the total number of coefficients used in the expansion. More...
 
int GetTotPoints () const
 This function returns the total number of quadrature points used in the element. More...
 
LibUtilities::BasisType GetBasisType (const int dir) const
 This function returns the type of basis used in the dir direction. More...
 
int GetBasisNumModes (const int dir) const
 This function returns the number of expansion modes in the dir direction. More...
 
int EvalBasisNumModesMax (void) const
 This function returns the maximum number of expansion modes over all local directions. More...
 
LibUtilities::PointsType GetPointsType (const int dir) const
 This function returns the type of quadrature points used in the dir direction. More...
 
int GetNumPoints (const int dir) const
 This function returns the number of quadrature points in the dir direction. More...
 
const Array< OneD, const NekDouble > & GetPoints (const int dir) const
 This function returns a pointer to the array containing the quadrature points in dir direction. More...
 
int GetNverts () const
 This function returns the number of vertices of the expansion domain. More...
 
int GetTraceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th trace. More...
 
int GetTraceIntNcoeffs (const int i) const
 
int GetTraceNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th trace. More...
 
const LibUtilities::BasisKey GetTraceBasisKey (const int i, int k=-1) const
 This function returns the basis key belonging to the i-th trace. More...
 
LibUtilities::PointsKey GetTracePointsKey (const int i, int k=-1) const
 This function returns the basis key belonging to the i-th trace. More...
 
int NumBndryCoeffs (void) const
 
int NumDGBndryCoeffs (void) const
 
const LibUtilities::PointsKey GetNodalPointsKey () const
 This function returns the type of expansion Nodal point type if defined. More...
 
int GetNtraces () 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...
 
std::shared_ptr< StdExpansionGetStdExp (void) const
 
std::shared_ptr< StdExpansionGetLinStdExp (void) const
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion ()
 
bool IsNodalNonTensorialExp ()
 
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space. More...
 
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space. More...
 
void FwdTrans_BndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble Integral (const Array< OneD, const NekDouble > &inarray)
 This function integrates the specified function over the domain. More...
 
void FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 This function fills the array outarray with the mode-th mode of the expansion. More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 this function calculates the inner product of a given function f with the different modes of the expansion More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check)
 
void IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase (const Array< OneD, const NekDouble > &direction, 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)
 
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)
 
int CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
NekDouble StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
int GetCoordim ()
 
void GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
void GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
int GetVertexMap (const int localVertexId, bool useCoeffPacking=false)
 
void GetTraceToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
 
void GetTraceInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eForwards)
 
void GetTraceNumModes (const int tid, int &numModes0, int &numModes1, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2)
 
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 ExponentialFilter (Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
 
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)
 
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...
 
NekDouble PhysEvaluateBasis (const Array< OneD, const NekDouble > &coords, int mode)
 This function evaluates the basis function mode mode at a point coords 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...
 
void LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi)
 Convert local collapsed coordinates eta into local cartesian coordinate xi. 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)
 
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 LibUtilities::PointsKeyVector GetPointsKeys () const
 
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 >
std::shared_ptr< T > as ()
 
void IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
void GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 
- Public Member Functions inherited from Nektar::LocalRegions::Expansion3D
 Expansion3D (SpatialDomains::Geometry3DSharedPtr pGeom)
 
virtual ~Expansion3D ()
 
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 v_ReOrientTracePhysMap (const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1)
 
void v_NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
Array< OneD, unsigned int > GetEdgeInverseBoundaryMap (int eid)
 
Array< OneD, unsigned int > GetTraceInverseBoundaryMap (int fid, StdRegions::Orientation faceOrient=StdRegions::eNoOrientation, int P1=-1, int P2=-1)
 
void GetInverseBoundaryMaps (Array< OneD, unsigned int > &vmap, Array< OneD, Array< OneD, unsigned int > > &emap, Array< OneD, Array< OneD, unsigned int > > &fmap)
 
- Public Member Functions inherited from Nektar::LocalRegions::Expansion
 Expansion (SpatialDomains::GeometrySharedPtr pGeom)
 
 Expansion (const Expansion &pSrc)
 
virtual ~Expansion ()
 
void SetTraceExp (const int traceid, ExpansionSharedPtr &f)
 
ExpansionSharedPtr GetTraceExp (const int traceid)
 
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 ()
 
IndexMapValuesSharedPtr CreateIndexMap (const IndexMapKey &ikey)
 
const SpatialDomains::GeomFactorsSharedPtrGetMetricInfo () const
 
DNekMatSharedPtr BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
 
DNekMatSharedPtr BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd)
 
void ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType)
 
void AddEdgeNormBoundaryInt (const int edge, const std::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 std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
void AddFaceNormBoundaryInt (const int face, const std::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)
 
NekDouble VectorFlux (const Array< OneD, Array< OneD, NekDouble > > &vec)
 
IndexMapValuesSharedPtr GetIndexMap (const IndexMapKey &ikey)
 
void AlignVectorToCollapsedDir (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
ExpansionSharedPtr GetLeftAdjacentElementExp () const
 
ExpansionSharedPtr GetRightAdjacentElementExp () const
 
int GetLeftAdjacentElementTrace () const
 
int GetRightAdjacentElementTrace () const
 
void SetAdjacentElementExp (int traceid, ExpansionSharedPtr &e)
 
StdRegions::Orientation GetTraceOrient (int trace)
 
void SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void DivideByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Divided by the metric jacobi and quadrature weights. More...
 
void GetTraceQFactors (const int trace, 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 GetTracePhysVals (const int trace, const StdRegions::StdExpansionSharedPtr &TraceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=StdRegions::eNoOrientation)
 
void GetTracePhysMap (const int edge, Array< OneD, int > &outarray)
 
void ReOrientTracePhysMap (const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1)
 
const NormalVectorGetTraceNormal (const int id) const
 
void ComputeTraceNormal (const int id)
 
const Array< OneD, const NekDouble > & GetPhysNormals (void)
 
void SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
void SetUpPhysNormals (const int trace)
 
void AddRobinMassMatrix (const int traceid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
virtual void AddRobinTraceContribution (const int traceid, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs)
 
const Array< OneD, const NekDouble > & GetElmtBndNormDirElmtLen (const int nbnd) const
 

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 void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Calculates the inner product \( I_{pqr} = (u, \partial_{x_i} \phi_{pqr}) \). More...
 
void v_IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_AlignVectorToCollapsedDir (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual StdRegions::StdExpansionSharedPtr v_GetStdExp (void) const
 
virtual StdRegions::StdExpansionSharedPtr v_GetLinStdExp (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 void v_ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int mode_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType)
 
NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
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_GetTracePhysMap (const int face, Array< OneD, int > &outarray)
 
void v_ComputeTraceNormal (const int face)
 
virtual void v_SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdRegions::StdMatrixKey &mkey)
 
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_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_LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 
virtual void v_LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi)
 
virtual void v_FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 
virtual void v_GetTraceNumModes (const int fid, int &numModes0, int &numModes1, Orientation faceOrient=eDir1FwdDir1_Dir2FwdDir2)
 
NekDouble v_PhysEvaluateBasis (const Array< OneD, const NekDouble > &coords, int mode) final
 
virtual int v_GetNverts () const
 
virtual int v_GetNedges () const
 
virtual int v_GetNtraces () const
 
virtual LibUtilities::ShapeType v_DetShapeType () const
 
virtual int v_NumBndryCoeffs () const
 
virtual int v_NumDGBndryCoeffs () const
 
virtual int v_GetTraceNcoeffs (const int i) const
 
virtual int v_GetTraceIntNcoeffs (const int i) const
 
virtual int v_GetTraceNumPoints (const int i) const
 
virtual int v_GetEdgeNcoeffs (const int i) const
 
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
virtual const LibUtilities::BasisKey v_GetTraceBasisKey (const int i, const int k) const
 
virtual int v_GetVertexMap (int localVertexId, bool useCoeffPacking=false)
 
virtual void v_GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
virtual void v_GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
virtual void v_GetTraceToElementMap (const int fid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation faceOrient, int P, int Q)
 
virtual void v_GetEdgeInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2)
 
virtual void v_GetTraceInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2)
 
virtual void v_ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion3D
virtual NekDouble v_PhysEvaluate (const Array< OneD, 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_GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 
- 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...
 
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 IProductWRTDirectionalDerivBase_SumFac (const Array< OneD, const NekDouble > &direction, 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 (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
 
template<int DIR>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals)
 This function performs the barycentric interpolation of the polynomial stored in coord at a point physvals using barycentric interpolation weights in direction. More...
 
template<int DIR>
NekDouble BaryEvaluateBasis (const NekDouble &coord, const int &mode)
 
- 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_GetTraceOrient (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)
 Extract the physical values along face face from inarray into outarray following the local face orientation and point distribution defined by defined in FaceExp. More...
 
void GetPhysFaceVarCoeffsFromElement (const int face, ExpansionSharedPtr &FaceExp, const Array< OneD, const NekDouble > &varcoeff, Array< OneD, NekDouble > &outarray)
 
virtual Array< OneD, NekDoublev_GetnFacecdotMF (const int dir, const int face, ExpansionSharedPtr &FaceExp_f, const Array< OneD, const Array< OneD, NekDouble > > &normals, const StdRegions::VarCoeffMap &varcoeffs)
 
virtual DNekMatSharedPtr v_BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
 
virtual DNekMatSharedPtr v_BuildInverseTransformationMatrix (const DNekScalMatSharedPtr &transformationmatrix)
 Build inverse and inverse transposed transformation matrix: \(\mathbf{R^{-1}}\) and \(\mathbf{R^{-T}}\). More...
 
virtual DNekMatSharedPtr v_BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd)
 
virtual const NormalVectorv_GetTraceNormal (const int face) const
 
- Protected Member Functions inherited from Nektar::LocalRegions::Expansion
void ComputeLaplacianMetric ()
 
void ComputeQuadratureMetric ()
 
void ComputeGmatcdotMF (const Array< TwoD, const NekDouble > &df, const Array< OneD, const NekDouble > &direction, Array< OneD, Array< OneD, NekDouble > > &dfdir)
 
virtual void v_MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_DivideByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
Array< OneD, NekDoublev_GetMF (const int dir, const int shapedim, const StdRegions::VarCoeffMap &varcoeffs)
 
Array< OneD, NekDoublev_GetMFDiv (const int dir, const StdRegions::VarCoeffMap &varcoeffs)
 
Array< OneD, NekDoublev_GetMFMag (const int dir, const StdRegions::VarCoeffMap &varcoeffs)
 
virtual void v_AddEdgeNormBoundaryInt (const int edge, const std::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 std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddFaceNormBoundaryInt (const int face, const std::shared_ptr< Expansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
virtual NekDouble v_VectorFlux (const Array< OneD, Array< OneD, NekDouble > > &vec)
 
virtual void v_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_GetTraceQFactors (const int trace, Array< OneD, NekDouble > &outarray)
 
virtual const Array< OneD, const NekDouble > & v_GetPhysNormals (void)
 
virtual void v_SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
virtual void v_SetUpPhysNormals (const int id)
 
virtual void v_AddRobinTraceContribution (const int traceid, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs)
 

Private Member Functions

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

Private Attributes

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

Additional Inherited Members

- Protected Attributes inherited from Nektar::StdRegions::StdExpansion
Array< OneD, LibUtilities::BasisSharedPtrm_base
 
int m_elmt_id
 
int m_ncoeffs
 
LibUtilities::NekManager< StdMatrixKey, DNekMat, StdMatrixKey::opLessm_stdMatrixManager
 
LibUtilities::NekManager< StdMatrixKey, DNekBlkMat, StdMatrixKey::opLessm_stdStaticCondMatrixManager
 
- Protected Attributes inherited from Nektar::LocalRegions::Expansion3D
std::map< int, NormalVectorm_faceNormals
 
- Protected Attributes inherited from Nektar::LocalRegions::Expansion
LibUtilities::NekManager< IndexMapKey, IndexMapValues, IndexMapKey::opLessm_indexMapManager
 
std::vector< ExpansionWeakPtrm_traceExp
 
SpatialDomains::GeometrySharedPtr m_geom
 
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
 
MetricMap m_metrics
 
ExpansionWeakPtr m_elementLeft
 
ExpansionWeakPtr m_elementRight
 
int m_elementTraceLeft = -1
 
int m_elementTraceRight = -1
 
std::map< int, Array< OneD, NekDouble > > m_elmtBndNormDirElmtLen
 the element length in each element boundary(Vertex, edge or face) normal direction calculated based on the local m_metricinfo times the standard element length (which is 2.0) More...
 

Detailed Description

Definition at line 50 of file PyrExp.h.

Constructor & Destructor Documentation

◆ PyrExp() [1/2]

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 45 of file PyrExp.cpp.

48  :
50  Ba.GetNumModes(),
51  Bb.GetNumModes(),
52  Bc.GetNumModes()),
53  3, Ba, Bb, Bc),
55  Ba.GetNumModes(),
56  Bb.GetNumModes(),
57  Bc.GetNumModes()),
58  Ba, Bb, Bc),
59  StdPyrExp (Ba,Bb,Bc),
60  Expansion (geom),
61  Expansion3D (geom),
63  std::bind(&PyrExp::CreateMatrix, this, std::placeholders::_1),
64  std::string("PyrExpMatrix")),
66  std::bind(&PyrExp::CreateStaticCondMatrix, this, std::placeholders::_1),
67  std::string("PyrExpStaticCondMatrix"))
68  {
69  }
Expansion3D(SpatialDomains::Geometry3DSharedPtr pGeom)
Definition: Expansion3D.h:61
Expansion(SpatialDomains::GeometrySharedPtr pGeom)
Definition: Expansion.cpp:47
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: PyrExp.h:186
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: PyrExp.h:187
DNekScalMatSharedPtr CreateMatrix(const MatrixKey &mkey)
Definition: PyrExp.cpp:1120
DNekScalBlkMatSharedPtr CreateStaticCondMatrix(const MatrixKey &mkey)
Definition: PyrExp.cpp:1268
StdExpansion()
Default Constructor.
int getNumberOfCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:240

◆ PyrExp() [2/2]

Nektar::LocalRegions::PyrExp::PyrExp ( const PyrExp T)

Definition at line 71 of file PyrExp.cpp.

71  :
72  StdExpansion (T),
73  StdExpansion3D(T),
74  StdPyrExp (T),
75  Expansion (T),
76  Expansion3D (T),
77  m_matrixManager(T.m_matrixManager),
78  m_staticCondMatrixManager(T.m_staticCondMatrixManager)
79  {
80  }

◆ ~PyrExp()

Nektar::LocalRegions::PyrExp::~PyrExp ( )

Definition at line 82 of file PyrExp.cpp.

83  {
84  }

Member Function Documentation

◆ CreateMatrix()

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

Definition at line 1120 of file PyrExp.cpp.

1121  {
1122  DNekScalMatSharedPtr returnval;
1124 
1125  ASSERTL2(m_metricinfo->GetGtype() != SpatialDomains::eNoGeomType,"Geometric information is not set up");
1126 
1127  switch(mkey.GetMatrixType())
1128  {
1129  case StdRegions::eMass:
1130  {
1131  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1132  {
1133  NekDouble one = 1.0;
1134  DNekMatSharedPtr mat = GenMatrix(mkey);
1135  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,mat);
1136  }
1137  else
1138  {
1139  NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
1140  DNekMatSharedPtr mat = GetStdMatrix(mkey);
1141  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(jac,mat);
1142  }
1143  }
1144  break;
1145  case StdRegions::eInvMass:
1146  {
1147  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1148  {
1149  NekDouble one = 1.0;
1150  StdRegions::StdMatrixKey masskey(StdRegions::eMass,DetShapeType(),
1151  *this);
1152  DNekMatSharedPtr mat = GenMatrix(masskey);
1153  mat->Invert();
1154  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,mat);
1155  }
1156  else
1157  {
1158  NekDouble fac = 1.0/(m_metricinfo->GetJac(ptsKeys))[0];
1159  DNekMatSharedPtr mat = GetStdMatrix(mkey);
1160  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(fac,mat);
1161  }
1162  }
1163  break;
1165  {
1166  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed ||
1167  mkey.GetNVarCoeff() > 0 ||
1168  mkey.ConstFactorExists(
1170  {
1171  NekDouble one = 1.0;
1172  DNekMatSharedPtr mat = GenMatrix(mkey);
1173 
1174  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,mat);
1175  }
1176  else
1177  {
1178  MatrixKey lap00key(StdRegions::eLaplacian00,
1179  mkey.GetShapeType(), *this);
1180  MatrixKey lap01key(StdRegions::eLaplacian01,
1181  mkey.GetShapeType(), *this);
1182  MatrixKey lap02key(StdRegions::eLaplacian02,
1183  mkey.GetShapeType(), *this);
1184  MatrixKey lap11key(StdRegions::eLaplacian11,
1185  mkey.GetShapeType(), *this);
1186  MatrixKey lap12key(StdRegions::eLaplacian12,
1187  mkey.GetShapeType(), *this);
1188  MatrixKey lap22key(StdRegions::eLaplacian22,
1189  mkey.GetShapeType(), *this);
1190 
1191  DNekMat &lap00 = *GetStdMatrix(lap00key);
1192  DNekMat &lap01 = *GetStdMatrix(lap01key);
1193  DNekMat &lap02 = *GetStdMatrix(lap02key);
1194  DNekMat &lap11 = *GetStdMatrix(lap11key);
1195  DNekMat &lap12 = *GetStdMatrix(lap12key);
1196  DNekMat &lap22 = *GetStdMatrix(lap22key);
1197 
1198  NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
1199  Array<TwoD, const NekDouble> gmat =
1200  m_metricinfo->GetGmat(ptsKeys);
1201 
1202  int rows = lap00.GetRows();
1203  int cols = lap00.GetColumns();
1204 
1206  ::AllocateSharedPtr(rows,cols);
1207 
1208  (*lap) = gmat[0][0]*lap00
1209  + gmat[4][0]*lap11
1210  + gmat[8][0]*lap22
1211  + gmat[3][0]*(lap01 + Transpose(lap01))
1212  + gmat[6][0]*(lap02 + Transpose(lap02))
1213  + gmat[7][0]*(lap12 + Transpose(lap12));
1214 
1215  returnval = MemoryManager<DNekScalMat>
1216  ::AllocateSharedPtr(jac, lap);
1217  }
1218  }
1219  break;
1221  {
1222  NekDouble factor = mkey.GetConstFactor(StdRegions::eFactorLambda);
1223  MatrixKey masskey(StdRegions::eMass, mkey.GetShapeType(), *this);
1224  DNekScalMat &MassMat = *(this->m_matrixManager[masskey]);
1225  MatrixKey lapkey(StdRegions::eLaplacian, mkey.GetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
1226  DNekScalMat &LapMat = *(this->m_matrixManager[lapkey]);
1227 
1228  int rows = LapMat.GetRows();
1229  int cols = LapMat.GetColumns();
1230 
1232 
1233  (*helm) = LapMat + factor*MassMat;
1234 
1235  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(1.0, helm);
1236  }
1237  break;
1239  {
1240  NekDouble one = 1.0;
1241  MatrixKey helmkey(StdRegions::eHelmholtz, mkey.GetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
1242  DNekScalBlkMatSharedPtr helmStatCond = GetLocStaticCondMatrix(helmkey);
1243  DNekScalMatSharedPtr A =helmStatCond->GetBlock(0,0);
1245 
1247  }
1248  break;
1250  {
1251  NekDouble one = 1.0;
1252  MatrixKey masskey(StdRegions::eMass, mkey.GetShapeType(), *this);
1253  DNekScalBlkMatSharedPtr massStatCond = GetLocStaticCondMatrix(masskey);
1254  DNekScalMatSharedPtr A =massStatCond->GetBlock(0,0);
1256 
1258  }
1259  break;
1260  default:
1261  NEKERROR(ErrorUtil::efatal, "Matrix creation not defined");
1262  break;
1263  }
1264 
1265  return returnval;
1266  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:274
DNekMatSharedPtr BuildVertexMatrix(const DNekScalMatSharedPtr &r_bnd)
Definition: Expansion.cpp:103
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:273
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:617
const LibUtilities::PointsKeyVector GetPointsKeys() const
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix(const LocalRegions::MatrixKey &mkey)
Definition: StdExpansion.h:660
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:376
DNekMatSharedPtr GenMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:850
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:246
@ eNoGeomType
No type defined.
@ eDeformed
Geometry is curved or has non-constant factors.
NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag > DNekScalMat
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
NekMatrix< NekDouble, StandardMatrixTag > DNekMat
Definition: NekTypeDefs.hpp:51
NekMatrix< InnerMatrixType, BlockMatrixTag > Transpose(NekMatrix< InnerMatrixType, BlockMatrixTag > &rhs)
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:73
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
double NekDouble

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL2, Nektar::LocalRegions::Expansion::BuildVertexMatrix(), Nektar::StdRegions::StdMatrixKey::ConstFactorExists(), Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::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::ePreconLinearSpace, Nektar::StdRegions::ePreconLinearSpaceMass, Nektar::StdRegions::StdExpansion::GenMatrix(), Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdMatrixKey::GetConstFactors(), Nektar::StdRegions::StdExpansion::GetLocStaticCondMatrix(), 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().

◆ CreateStaticCondMatrix()

DNekScalBlkMatSharedPtr Nektar::LocalRegions::PyrExp::CreateStaticCondMatrix ( const MatrixKey mkey)
protected

Definition at line 1268 of file PyrExp.cpp.

1269  {
1270  DNekScalBlkMatSharedPtr returnval;
1271 
1272  ASSERTL2(m_metricinfo->GetGtype() != SpatialDomains::eNoGeomType,"Geometric information is not set up");
1273 
1274  // set up block matrix system
1275  unsigned int nbdry = NumBndryCoeffs();
1276  unsigned int nint = (unsigned int)(m_ncoeffs - nbdry);
1277  unsigned int exp_size[] = {nbdry, nint};
1278  unsigned int nblks = 2;
1279  returnval = MemoryManager<DNekScalBlkMat>::AllocateSharedPtr(nblks, nblks, exp_size, exp_size); //Really need a constructor which takes Arrays
1280  NekDouble factor = 1.0;
1281 
1282  switch(mkey.GetMatrixType())
1283  {
1285  case StdRegions::eHelmholtz: // special case since Helmholtz not defined in StdRegions
1286 
1287  // use Deformed case for both regular and deformed geometries
1288  factor = 1.0;
1289  goto UseLocRegionsMatrix;
1290  break;
1291  default:
1292  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1293  {
1294  factor = 1.0;
1295  goto UseLocRegionsMatrix;
1296  }
1297  else
1298  {
1299  DNekScalMatSharedPtr mat = GetLocMatrix(mkey);
1300  factor = mat->Scale();
1301  goto UseStdRegionsMatrix;
1302  }
1303  break;
1304  UseStdRegionsMatrix:
1305  {
1306  NekDouble invfactor = 1.0/factor;
1307  NekDouble one = 1.0;
1309  DNekScalMatSharedPtr Atmp;
1310  DNekMatSharedPtr Asubmat;
1311 
1312  //TODO: check below
1313  returnval->SetBlock(0,0,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(factor,Asubmat = mat->GetBlock(0,0)));
1314  returnval->SetBlock(0,1,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,Asubmat = mat->GetBlock(0,1)));
1315  returnval->SetBlock(1,0,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(factor,Asubmat = mat->GetBlock(1,0)));
1316  returnval->SetBlock(1,1,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(invfactor,Asubmat = mat->GetBlock(1,1)));
1317  }
1318  break;
1319  UseLocRegionsMatrix:
1320  {
1321  int i,j;
1322  NekDouble invfactor = 1.0/factor;
1323  NekDouble one = 1.0;
1324  DNekScalMat &mat = *GetLocMatrix(mkey);
1329 
1330  Array<OneD,unsigned int> bmap(nbdry);
1331  Array<OneD,unsigned int> imap(nint);
1332  GetBoundaryMap(bmap);
1333  GetInteriorMap(imap);
1334 
1335  for(i = 0; i < nbdry; ++i)
1336  {
1337  for(j = 0; j < nbdry; ++j)
1338  {
1339  (*A)(i,j) = mat(bmap[i],bmap[j]);
1340  }
1341 
1342  for(j = 0; j < nint; ++j)
1343  {
1344  (*B)(i,j) = mat(bmap[i],imap[j]);
1345  }
1346  }
1347 
1348  for(i = 0; i < nint; ++i)
1349  {
1350  for(j = 0; j < nbdry; ++j)
1351  {
1352  (*C)(i,j) = mat(imap[i],bmap[j]);
1353  }
1354 
1355  for(j = 0; j < nint; ++j)
1356  {
1357  (*D)(i,j) = mat(imap[i],imap[j]);
1358  }
1359  }
1360 
1361  // Calculate static condensed system
1362  if(nint)
1363  {
1364  D->Invert();
1365  (*B) = (*B)*(*D);
1366  (*A) = (*A) - (*B)*(*C);
1367  }
1368 
1369  DNekScalMatSharedPtr Atmp;
1370 
1371  returnval->SetBlock(0,0,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(factor,A));
1372  returnval->SetBlock(0,1,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,B));
1373  returnval->SetBlock(1,0,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(factor,C));
1374  returnval->SetBlock(1,1,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(invfactor,D));
1375  }
1376  }
1377  return returnval;
1378  }
DNekScalMatSharedPtr GetLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:90
void GetBoundaryMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:687
DNekBlkMatSharedPtr GetStdStaticCondMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:622
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:692
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:71

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

◆ v_AlignVectorToCollapsedDir()

void Nektar::LocalRegions::PyrExp::v_AlignVectorToCollapsedDir ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 410 of file PyrExp.cpp.

414  {
415  const int nquad0 = m_base[0]->GetNumPoints();
416  const int nquad1 = m_base[1]->GetNumPoints();
417  const int nquad2 = m_base[2]->GetNumPoints();
418  const int order0 = m_base[0]->GetNumModes ();
419  const int order1 = m_base[1]->GetNumModes ();
420  const int nqtot = nquad0*nquad1*nquad2;
421 
422  const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
423  const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
424  const Array<OneD, const NekDouble> &z2 = m_base[2]->GetZ();
425 
426  Array<OneD, NekDouble> gfac0(nquad0 );
427  Array<OneD, NekDouble> gfac1(nquad1 );
428  Array<OneD, NekDouble> gfac2(nquad2 );
429  Array<OneD, NekDouble> tmp5 (nqtot );
430  Array<OneD, NekDouble> wsp (std::max(nqtot,
431  order0*nquad2*(nquad1+order1)));
432 
433  Array<OneD, NekDouble> tmp2 = outarray[0];
434  Array<OneD, NekDouble> tmp3 = outarray[1];
435  Array<OneD, NekDouble> tmp4 = outarray[2];
436 
437  const Array<TwoD, const NekDouble>& df =
438  m_metricinfo->GetDerivFactors(GetPointsKeys());
439 
440  Array<OneD, NekDouble> tmp1;
441  tmp1 = inarray;
442 
443  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
444  {
445  Vmath::Vmul(nqtot,&df[3*dir][0], 1,tmp1.get(),1,tmp2.get(),1);
446  Vmath::Vmul(nqtot,&df[3*dir+1][0],1,tmp1.get(),1,tmp3.get(),1);
447  Vmath::Vmul(nqtot,&df[3*dir+2][0],1,tmp1.get(),1,tmp4.get(),1);
448  }
449  else
450  {
451  Vmath::Smul(nqtot, df[3*dir][0], tmp1.get(),1,tmp2.get(), 1);
452  Vmath::Smul(nqtot, df[3*dir+1][0],tmp1.get(),1,tmp3.get(), 1);
453  Vmath::Smul(nqtot, df[3*dir+2][0],tmp1.get(),1,tmp4.get(), 1);
454  }
455 
456  // set up geometric factor: (1+z0)/2
457  for (int i = 0; i < nquad0; ++i)
458  {
459  gfac0[i] = 0.5*(1+z0[i]);
460  }
461 
462  // set up geometric factor: (1+z1)/2
463  for(int i = 0; i < nquad1; ++i)
464  {
465  gfac1[i] = 0.5*(1+z1[i]);
466  }
467 
468  // Set up geometric factor: 2/(1-z2)
469  for (int i = 0; i < nquad2; ++i)
470  {
471  gfac2[i] = 2.0/(1-z2[i]);
472  }
473 
474  const int nq01 = nquad0*nquad1;
475 
476  for (int i = 0; i < nquad2; ++i)
477  {
478  Vmath::Smul(nq01,gfac2[i],&tmp2[0]+i*nq01,1,&tmp2[0]+i*nq01,1); // 2/(1-z2) for d/dxi_0
479  Vmath::Smul(nq01,gfac2[i],&tmp3[0]+i*nq01,1,&tmp3[0]+i*nq01,1); // 2/(1-z2) for d/dxi_1
480  Vmath::Smul(nq01,gfac2[i],&tmp4[0]+i*nq01,1,&tmp5[0]+i*nq01,1); // 2/(1-z2) for d/dxi_2
481  }
482 
483  // (1+z0)/(1-z2) for d/d eta_0
484  for (int i = 0; i < nquad1*nquad2; ++i)
485  {
486  Vmath::Vmul(nquad0,&gfac0[0],1,
487  &tmp5[0]+i*nquad0,1,
488  &wsp[0]+i*nquad0,1);
489  }
490 
491  Vmath::Vadd(nqtot, &tmp2[0], 1, &wsp[0], 1, &tmp2[0], 1);
492 
493  // (1+z1)/(1-z2) for d/d eta_1
494  for (int i = 0; i < nquad1*nquad2; ++i)
495  {
496  Vmath::Smul(nquad0,gfac1[i%nquad1],
497  &tmp5[0]+i*nquad0,1,
498  &tmp5[0]+i*nquad0,1);
499  }
500  Vmath::Vadd(nqtot, &tmp3[0], 1, &tmp5[0], 1, &tmp3[0], 1);
501  }
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:192
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:322
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:225

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

Referenced by v_IProductWRTDerivBase_SumFac().

◆ v_ComputeLaplacianMetric()

void Nektar::LocalRegions::PyrExp::v_ComputeLaplacianMetric ( )
protectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1380 of file PyrExp.cpp.

1381  {
1382  if (m_metrics.count(eMetricQuadrature) == 0)
1383  {
1385  }
1386 
1387  int i, j;
1388  const unsigned int nqtot = GetTotPoints();
1389  const unsigned int dim = 3;
1390  const MetricType m[3][3] = {
1394  };
1395 
1396  for (unsigned int i = 0; i < dim; ++i)
1397  {
1398  for (unsigned int j = i; j < dim; ++j)
1399  {
1400  m_metrics[m[i][j]] = Array<OneD, NekDouble>(nqtot);
1401  }
1402  }
1403 
1404  // Define shorthand synonyms for m_metrics storage
1405  Array<OneD,NekDouble> g0 (m_metrics[m[0][0]]);
1406  Array<OneD,NekDouble> g1 (m_metrics[m[1][1]]);
1407  Array<OneD,NekDouble> g2 (m_metrics[m[2][2]]);
1408  Array<OneD,NekDouble> g3 (m_metrics[m[0][1]]);
1409  Array<OneD,NekDouble> g4 (m_metrics[m[0][2]]);
1410  Array<OneD,NekDouble> g5 (m_metrics[m[1][2]]);
1411 
1412  // Allocate temporary storage
1413  Array<OneD,NekDouble> alloc(9*nqtot,0.0);
1414  Array<OneD,NekDouble> h0 (nqtot, alloc);
1415  Array<OneD,NekDouble> h1 (nqtot, alloc+ 1*nqtot);
1416  Array<OneD,NekDouble> h2 (nqtot, alloc+ 2*nqtot);
1417  Array<OneD,NekDouble> wsp1 (nqtot, alloc+ 3*nqtot);
1418  Array<OneD,NekDouble> wsp2 (nqtot, alloc+ 4*nqtot);
1419  Array<OneD,NekDouble> wsp3 (nqtot, alloc+ 5*nqtot);
1420  Array<OneD,NekDouble> wsp4 (nqtot, alloc+ 6*nqtot);
1421  Array<OneD,NekDouble> wsp5 (nqtot, alloc+ 7*nqtot);
1422  Array<OneD,NekDouble> wsp6 (nqtot, alloc+ 8*nqtot);
1423 
1424  const Array<TwoD, const NekDouble>& df =
1425  m_metricinfo->GetDerivFactors(GetPointsKeys());
1426  const Array<OneD, const NekDouble>& z0 = m_base[0]->GetZ();
1427  const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
1428  const Array<OneD, const NekDouble>& z2 = m_base[2]->GetZ();
1429  const unsigned int nquad0 = m_base[0]->GetNumPoints();
1430  const unsigned int nquad1 = m_base[1]->GetNumPoints();
1431  const unsigned int nquad2 = m_base[2]->GetNumPoints();
1432 
1433  // Populate collapsed coordinate arrays h0, h1 and h2.
1434  for(j = 0; j < nquad2; ++j)
1435  {
1436  for(i = 0; i < nquad1; ++i)
1437  {
1438  Vmath::Fill(nquad0, 2.0/(1.0-z2[j]), &h0[0]+i*nquad0 + j*nquad0*nquad1,1);
1439  Vmath::Fill(nquad0, 1.0/(1.0-z2[j]), &h1[0]+i*nquad0 + j*nquad0*nquad1,1);
1440  Vmath::Fill(nquad0, (1.0+z1[i])/(1.0-z2[j]), &h2[0]+i*nquad0 + j*nquad0*nquad1,1);
1441  }
1442  }
1443  for(i = 0; i < nquad0; i++)
1444  {
1445  Blas::Dscal(nquad1*nquad2, 1+z0[i], &h1[0]+i, nquad0);
1446  }
1447 
1448  // Step 3. Construct combined metric terms for physical space to
1449  // collapsed coordinate system.
1450  // Order of construction optimised to minimise temporary storage
1451  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1452  {
1453  // f_{1k}
1454  Vmath::Vvtvvtp(nqtot, &df[0][0], 1, &h0[0], 1, &df[2][0], 1, &h1[0], 1, &wsp1[0], 1);
1455  Vmath::Vvtvvtp(nqtot, &df[3][0], 1, &h0[0], 1, &df[5][0], 1, &h1[0], 1, &wsp2[0], 1);
1456  Vmath::Vvtvvtp(nqtot, &df[6][0], 1, &h0[0], 1, &df[8][0], 1, &h1[0], 1, &wsp3[0], 1);
1457 
1458  // g0
1459  Vmath::Vvtvvtp(nqtot, &wsp1[0], 1, &wsp1[0], 1, &wsp2[0], 1, &wsp2[0], 1, &g0[0], 1);
1460  Vmath::Vvtvp (nqtot, &wsp3[0], 1, &wsp3[0], 1, &g0[0], 1, &g0[0], 1);
1461 
1462  // g4
1463  Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &wsp1[0], 1, &df[5][0], 1, &wsp2[0], 1, &g4[0], 1);
1464  Vmath::Vvtvp (nqtot, &df[8][0], 1, &wsp3[0], 1, &g4[0], 1, &g4[0], 1);
1465 
1466  // f_{2k}
1467  Vmath::Vvtvvtp(nqtot, &df[1][0], 1, &h0[0], 1, &df[2][0], 1, &h2[0], 1, &wsp4[0], 1);
1468  Vmath::Vvtvvtp(nqtot, &df[4][0], 1, &h0[0], 1, &df[5][0], 1, &h2[0], 1, &wsp5[0], 1);
1469  Vmath::Vvtvvtp(nqtot, &df[7][0], 1, &h0[0], 1, &df[8][0], 1, &h2[0], 1, &wsp6[0], 1);
1470 
1471  // g1
1472  Vmath::Vvtvvtp(nqtot, &wsp4[0], 1, &wsp4[0], 1, &wsp5[0], 1, &wsp5[0], 1, &g1[0], 1);
1473  Vmath::Vvtvp (nqtot, &wsp6[0], 1, &wsp6[0], 1, &g1[0], 1, &g1[0], 1);
1474 
1475  // g3
1476  Vmath::Vvtvvtp(nqtot, &wsp1[0], 1, &wsp4[0], 1, &wsp2[0], 1, &wsp5[0], 1, &g3[0], 1);
1477  Vmath::Vvtvp (nqtot, &wsp3[0], 1, &wsp6[0], 1, &g3[0], 1, &g3[0], 1);
1478 
1479  // g5
1480  Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &wsp4[0], 1, &df[5][0], 1, &wsp5[0], 1, &g5[0], 1);
1481  Vmath::Vvtvp (nqtot, &df[8][0], 1, &wsp6[0], 1, &g5[0], 1, &g5[0], 1);
1482 
1483  // g2
1484  Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &df[2][0], 1, &df[5][0], 1, &df[5][0], 1, &g2[0], 1);
1485  Vmath::Vvtvp (nqtot, &df[8][0], 1, &df[8][0], 1, &g2[0], 1, &g2[0], 1);
1486  }
1487  else
1488  {
1489  // f_{1k}
1490  Vmath::Svtsvtp(nqtot, df[0][0], &h0[0], 1, df[2][0], &h1[0], 1, &wsp1[0], 1);
1491  Vmath::Svtsvtp(nqtot, df[3][0], &h0[0], 1, df[5][0], &h1[0], 1, &wsp2[0], 1);
1492  Vmath::Svtsvtp(nqtot, df[6][0], &h0[0], 1, df[8][0], &h1[0], 1, &wsp3[0], 1);
1493 
1494  // g0
1495  Vmath::Vvtvvtp(nqtot, &wsp1[0], 1, &wsp1[0], 1, &wsp2[0], 1, &wsp2[0], 1, &g0[0], 1);
1496  Vmath::Vvtvp (nqtot, &wsp3[0], 1, &wsp3[0], 1, &g0[0], 1, &g0[0], 1);
1497 
1498  // g4
1499  Vmath::Svtsvtp(nqtot, df[2][0], &wsp1[0], 1, df[5][0], &wsp2[0], 1, &g4[0], 1);
1500  Vmath::Svtvp (nqtot, df[8][0], &wsp3[0], 1, &g4[0], 1, &g4[0], 1);
1501 
1502  // f_{2k}
1503  Vmath::Svtsvtp(nqtot, df[1][0], &h0[0], 1, df[2][0], &h2[0], 1, &wsp4[0], 1);
1504  Vmath::Svtsvtp(nqtot, df[4][0], &h0[0], 1, df[5][0], &h2[0], 1, &wsp5[0], 1);
1505  Vmath::Svtsvtp(nqtot, df[7][0], &h0[0], 1, df[8][0], &h2[0], 1, &wsp6[0], 1);
1506 
1507  // g1
1508  Vmath::Vvtvvtp(nqtot, &wsp4[0], 1, &wsp4[0], 1, &wsp5[0], 1, &wsp5[0], 1, &g1[0], 1);
1509  Vmath::Vvtvp (nqtot, &wsp6[0], 1, &wsp6[0], 1, &g1[0], 1, &g1[0], 1);
1510 
1511  // g3
1512  Vmath::Vvtvvtp(nqtot, &wsp1[0], 1, &wsp4[0], 1, &wsp2[0], 1, &wsp5[0], 1, &g3[0], 1);
1513  Vmath::Vvtvp (nqtot, &wsp3[0], 1, &wsp6[0], 1, &g3[0], 1, &g3[0], 1);
1514 
1515  // g5
1516  Vmath::Svtsvtp(nqtot, df[2][0], &wsp4[0], 1, df[5][0], &wsp5[0], 1, &g5[0], 1);
1517  Vmath::Svtvp (nqtot, df[8][0], &wsp6[0], 1, &g5[0], 1, &g5[0], 1);
1518 
1519  // g2
1520  Vmath::Fill(nqtot, df[2][0]*df[2][0] + df[5][0]*df[5][0] + df[8][0]*df[8][0], &g2[0], 1);
1521  }
1522 
1523  for (unsigned int i = 0; i < dim; ++i)
1524  {
1525  for (unsigned int j = i; j < dim; ++j)
1526  {
1528  m_metrics[m[i][j]]);
1529 
1530  }
1531  }
1532  }
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:134
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:733
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition: Blas.hpp:182
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:691
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:565
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:513
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45
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:625

References Nektar::LocalRegions::Expansion::ComputeQuadratureMetric(), Blas::Dscal(), 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().

◆ v_ComputeTraceNormal()

void Nektar::LocalRegions::PyrExp::v_ComputeTraceNormal ( const int  face)
protectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 750 of file PyrExp.cpp.

751  {
752  const SpatialDomains::GeomFactorsSharedPtr &geomFactors =
753  GetGeom()->GetMetricInfo();
754 
756  for(int i = 0; i < ptsKeys.size(); ++i)
757  {
758  // Need at least 2 points for computing normals
759  if (ptsKeys[i].GetNumPoints() == 1)
760  {
761  LibUtilities::PointsKey pKey(2, ptsKeys[i].GetPointsType());
762  ptsKeys[i] = pKey;
763  }
764  }
765 
766  SpatialDomains::GeomType type = geomFactors->GetGtype();
767  const Array<TwoD, const NekDouble> &df = geomFactors->GetDerivFactors(ptsKeys);
768  const Array<OneD, const NekDouble> &jac = geomFactors->GetJac(ptsKeys);
769 
770  LibUtilities::BasisKey tobasis0 = GetTraceBasisKey(face,0);
771  LibUtilities::BasisKey tobasis1 = GetTraceBasisKey(face,1);
772 
773  // Number of quadrature points in face expansion.
774  int nq_face = tobasis0.GetNumPoints()*tobasis1.GetNumPoints();
775 
776  int vCoordDim = GetCoordim();
777  int i;
778 
779  m_faceNormals[face] = Array<OneD, Array<OneD, NekDouble> >(vCoordDim);
780  Array<OneD, Array<OneD, NekDouble> > &normal = m_faceNormals[face];
781  for (i = 0; i < vCoordDim; ++i)
782  {
783  normal[i] = Array<OneD, NekDouble>(nq_face);
784  }
785 
786  size_t nqb = nq_face;
787  size_t nbnd= face;
788  m_elmtBndNormDirElmtLen[nbnd] = Array<OneD, NekDouble> {nqb, 0.0};
789  Array<OneD, NekDouble> &length = m_elmtBndNormDirElmtLen[nbnd];
790 
791  // Regular geometry case
792  if (type == SpatialDomains::eRegular ||
794  {
795  NekDouble fac;
796  // Set up normals
797  switch(face)
798  {
799  case 0:
800  {
801  for(i = 0; i < vCoordDim; ++i)
802  {
803  normal[i][0] = -df[3*i+2][0];
804  }
805  break;
806  }
807  case 1:
808  {
809  for(i = 0; i < vCoordDim; ++i)
810  {
811  normal[i][0] = -df[3*i+1][0];
812  }
813  break;
814  }
815  case 2:
816  {
817  for(i = 0; i < vCoordDim; ++i)
818  {
819  normal[i][0] = df[3*i][0]+df[3*i+2][0];
820  }
821  break;
822  }
823  case 3:
824  {
825  for(i = 0; i < vCoordDim; ++i)
826  {
827  normal[i][0] = df[3*i+1][0]+df[3*i+2][0];
828  }
829  break;
830  }
831  case 4:
832  {
833  for(i = 0; i < vCoordDim; ++i)
834  {
835  normal[i][0] = -df[3*i][0];
836  }
837  break;
838  }
839  default:
840  ASSERTL0(false,"face is out of range (face < 4)");
841  }
842 
843  // Normalise resulting vector.
844  fac = 0.0;
845  for(i = 0; i < vCoordDim; ++i)
846  {
847  fac += normal[i][0]*normal[i][0];
848  }
849  fac = 1.0/sqrt(fac);
850 
851  Vmath::Fill(nqb, fac, length, 1);
852 
853  for (i = 0; i < vCoordDim; ++i)
854  {
855  Vmath::Fill(nq_face,fac*normal[i][0],normal[i],1);
856  }
857  }
858  else
859  {
860  // Set up deformed normals.
861  int j, k;
862 
863  int nq0 = ptsKeys[0].GetNumPoints();
864  int nq1 = ptsKeys[1].GetNumPoints();
865  int nq2 = ptsKeys[2].GetNumPoints();
866  int nq01 = nq0*nq1;
867  int nqtot;
868 
869  // Determine number of quadrature points on the face.
870  if (face == 0)
871  {
872  nqtot = nq0*nq1;
873  }
874  else if (face == 1 || face == 3)
875  {
876  nqtot = nq0*nq2;
877  }
878  else
879  {
880  nqtot = nq1*nq2;
881  }
882 
883  LibUtilities::PointsKey points0;
884  LibUtilities::PointsKey points1;
885 
886  Array<OneD, NekDouble> faceJac(nqtot);
887  Array<OneD, NekDouble> normals(vCoordDim*nqtot,0.0);
888 
889  // Extract Jacobian along face and recover local derivatives
890  // (dx/dr) for polynomial interpolation by multiplying m_gmat by
891  // jacobian
892  switch(face)
893  {
894  case 0:
895  {
896  for(j = 0; j < nq01; ++j)
897  {
898  normals[j] = -df[2][j]*jac[j];
899  normals[nqtot+j] = -df[5][j]*jac[j];
900  normals[2*nqtot+j] = -df[8][j]*jac[j];
901  faceJac[j] = jac[j];
902  }
903 
904  points0 = ptsKeys[0];
905  points1 = ptsKeys[1];
906  break;
907  }
908 
909  case 1:
910  {
911  for (j = 0; j < nq0; ++j)
912  {
913  for(k = 0; k < nq2; ++k)
914  {
915  int tmp = j+nq01*k;
916  normals[j+k*nq0] =
917  -df[1][tmp]*jac[tmp];
918  normals[nqtot+j+k*nq0] =
919  -df[4][tmp]*jac[tmp];
920  normals[2*nqtot+j+k*nq0] =
921  -df[7][tmp]*jac[tmp];
922  faceJac[j+k*nq0] = jac[tmp];
923  }
924  }
925 
926  points0 = ptsKeys[0];
927  points1 = ptsKeys[2];
928  break;
929  }
930 
931  case 2:
932  {
933  for (j = 0; j < nq1; ++j)
934  {
935  for(k = 0; k < nq2; ++k)
936  {
937  int tmp = nq0-1+nq0*j+nq01*k;
938  normals[j+k*nq1] =
939  (df[0][tmp]+df[2][tmp])*jac[tmp];
940  normals[nqtot+j+k*nq1] =
941  (df[3][tmp]+df[5][tmp])*jac[tmp];
942  normals[2*nqtot+j+k*nq1] =
943  (df[6][tmp]+df[8][tmp])*jac[tmp];
944  faceJac[j+k*nq1] = jac[tmp];
945  }
946  }
947 
948  points0 = ptsKeys[1];
949  points1 = ptsKeys[2];
950  break;
951  }
952 
953  case 3:
954  {
955  for (j = 0; j < nq0; ++j)
956  {
957  for(k = 0; k < nq2; ++k)
958  {
959  int tmp = nq0*(nq1-1) + j + nq01*k;
960  normals[j+k*nq0] =
961  (df[1][tmp]+df[2][tmp])*jac[tmp];
962  normals[nqtot+j+k*nq0] =
963  (df[4][tmp]+df[5][tmp])*jac[tmp];
964  normals[2*nqtot+j+k*nq0] =
965  (df[7][tmp]+df[8][tmp])*jac[tmp];
966  faceJac[j+k*nq0] = jac[tmp];
967  }
968  }
969 
970  points0 = ptsKeys[0];
971  points1 = ptsKeys[2];
972  break;
973  }
974 
975  case 4:
976  {
977  for (j = 0; j < nq1; ++j)
978  {
979  for(k = 0; k < nq2; ++k)
980  {
981  int tmp = j*nq0+nq01*k;
982  normals[j+k*nq1] =
983  -df[0][tmp]*jac[tmp];
984  normals[nqtot+j+k*nq1] =
985  -df[3][tmp]*jac[tmp];
986  normals[2*nqtot+j+k*nq1] =
987  -df[6][tmp]*jac[tmp];
988  faceJac[j+k*nq1] = jac[tmp];
989  }
990  }
991 
992  points0 = ptsKeys[1];
993  points1 = ptsKeys[2];
994  break;
995  }
996 
997  default:
998  ASSERTL0(false,"face is out of range (face < 4)");
999  }
1000 
1001  Array<OneD, NekDouble> work (nq_face, 0.0);
1002  // Interpolate Jacobian and invert
1003  LibUtilities::Interp2D(points0, points1, faceJac,
1004  tobasis0.GetPointsKey(),
1005  tobasis1.GetPointsKey(),
1006  work);
1007  Vmath::Sdiv(nq_face, 1.0, &work[0], 1, &work[0], 1);
1008 
1009  // Interpolate normal and multiply by inverse Jacobian.
1010  for(i = 0; i < vCoordDim; ++i)
1011  {
1012  LibUtilities::Interp2D(points0, points1,
1013  &normals[i*nqtot],
1014  tobasis0.GetPointsKey(),
1015  tobasis1.GetPointsKey(),
1016  &normal[i][0]);
1017  Vmath::Vmul(nq_face,work,1,normal[i],1,normal[i],1);
1018  }
1019 
1020  // Normalise to obtain unit normals.
1021  Vmath::Zero(nq_face,work,1);
1022  for(i = 0; i < GetCoordim(); ++i)
1023  {
1024  Vmath::Vvtvp(nq_face,normal[i],1,normal[i],1,work,1,work,1);
1025  }
1026 
1027  Vmath::Vsqrt(nq_face,work,1,work,1);
1028  Vmath::Sdiv (nq_face,1.0,work,1,work,1);
1029 
1030  Vmath::Vcopy(nqb, work, 1, length, 1);
1031 
1032  for(i = 0; i < GetCoordim(); ++i)
1033  {
1034  Vmath::Vmul(nq_face,normal[i],1,work,1,normal[i],1);
1035  }
1036  }
1037  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
std::map< int, NormalVector > m_faceNormals
Definition: Expansion3D.h:123
std::map< int, Array< OneD, NekDouble > > m_elmtBndNormDirElmtLen
the element length in each element boundary(Vertex, edge or face) normal direction calculated based o...
Definition: Expansion.h:284
SpatialDomains::GeometrySharedPtr GetGeom() const
Definition: Expansion.cpp:172
const LibUtilities::BasisKey GetTraceBasisKey(const int i, int k=-1) const
This function returns the basis key belonging to the i-th trace.
Definition: StdExpansion.h:304
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:208
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:221
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:115
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition: GeomFactors.h:62
GeomType
Indicates the type of element geometry.
@ eRegular
Geometry is straight-sided with constant geometric factors.
@ eMovingRegular
Currently unused.
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
Definition: Vmath.cpp:475
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:291
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:436
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:267

References ASSERTL0, Nektar::SpatialDomains::eMovingRegular, Nektar::SpatialDomains::eRegular, Vmath::Fill(), Nektar::StdRegions::StdExpansion::GetCoordim(), Nektar::LocalRegions::Expansion::GetGeom(), Nektar::LibUtilities::BasisKey::GetNumPoints(), Nektar::StdRegions::StdExpansion::GetNumPoints(), Nektar::LibUtilities::BasisKey::GetPointsKey(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::StdRegions::StdExpansion::GetTraceBasisKey(), Nektar::LibUtilities::Interp2D(), Nektar::LocalRegions::Expansion::m_elmtBndNormDirElmtLen, Nektar::LocalRegions::Expansion3D::m_faceNormals, Vmath::Sdiv(), tinysimd::sqrt(), Vmath::Vcopy(), Vmath::Vmul(), Vmath::Vsqrt(), Vmath::Vvtvp(), and Vmath::Zero().

◆ v_CreateStdMatrix()

DNekMatSharedPtr Nektar::LocalRegions::PyrExp::v_CreateStdMatrix ( const StdRegions::StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdPyrExp.

Definition at line 1094 of file PyrExp.cpp.

1095  {
1096  LibUtilities::BasisKey bkey0 = m_base[0]->GetBasisKey();
1097  LibUtilities::BasisKey bkey1 = m_base[1]->GetBasisKey();
1098  LibUtilities::BasisKey bkey2 = m_base[2]->GetBasisKey();
1100  MemoryManager<StdPyrExp>::AllocateSharedPtr(bkey0, bkey1, bkey2);
1101 
1102  return tmp->GetStdMatrix(mkey);
1103  }
std::shared_ptr< StdPyrExp > StdPyrExpSharedPtr
Definition: StdPyrExp.h:278

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

◆ v_DropLocStaticCondMatrix()

void Nektar::LocalRegions::PyrExp::v_DropLocStaticCondMatrix ( const MatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1115 of file PyrExp.cpp.

1116  {
1117  m_staticCondMatrixManager.DeleteObject(mkey);
1118  }

References m_staticCondMatrixManager.

◆ v_ExtractDataToCoeffs()

void Nektar::LocalRegions::PyrExp::v_ExtractDataToCoeffs ( const NekDouble data,
const std::vector< unsigned int > &  nummodes,
const int  mode_offset,
NekDouble coeffs,
std::vector< LibUtilities::BasisType > &  fromType 
)
protectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 559 of file PyrExp.cpp.

565  {
566  int data_order0 = nummodes[mode_offset];
567  int fillorder0 = min(m_base[0]->GetNumModes(),data_order0);
568  int data_order1 = nummodes[mode_offset+1];
569  int order1 = m_base[1]->GetNumModes();
570  int fillorder1 = min(order1,data_order1);
571  int data_order2 = nummodes[mode_offset+2];
572  int order2 = m_base[2]->GetNumModes();
573  int fillorder2 = min(order2,data_order2);
574 
575  // Check if not same order or basis and if not make temp
576  // element to read in data
577  if (fromType[0] != m_base[0]->GetBasisType() ||
578  fromType[1] != m_base[1]->GetBasisType() ||
579  fromType[2] != m_base[2]->GetBasisType() ||
580  data_order0 != fillorder0 ||
581  data_order1 != fillorder1 ||
582  data_order2 != fillorder2)
583  {
584  // Construct a pyr with the appropriate basis type at our
585  // quadrature points, and one more to do a forwards
586  // transform. We can then copy the output to coeffs.
587  StdRegions::StdPyrExp tmpPyr(
588  LibUtilities::BasisKey(
589  fromType[0], data_order0, m_base[0]->GetPointsKey()),
590  LibUtilities::BasisKey(
591  fromType[1], data_order1, m_base[1]->GetPointsKey()),
592  LibUtilities::BasisKey(
593  fromType[2], data_order2, m_base[2]->GetPointsKey()));
594 
595  StdRegions::StdPyrExp tmpPyr2(m_base[0]->GetBasisKey(),
596  m_base[1]->GetBasisKey(),
597  m_base[2]->GetBasisKey());
598 
599  Array<OneD, const NekDouble> tmpData(tmpPyr.GetNcoeffs(), data);
600  Array<OneD, NekDouble> tmpBwd(tmpPyr2.GetTotPoints());
601  Array<OneD, NekDouble> tmpOut(tmpPyr2.GetNcoeffs());
602 
603  tmpPyr.BwdTrans(tmpData, tmpBwd);
604  tmpPyr2.FwdTrans(tmpBwd, tmpOut);
605  Vmath::Vcopy(tmpOut.size(), &tmpOut[0], 1, coeffs, 1);
606 
607  }
608  else
609  {
610  Vmath::Vcopy(m_ncoeffs,&data[0],1,coeffs,1);
611  }
612  }
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:158

References Nektar::StdRegions::StdExpansion::BwdTrans(), Nektar::StdRegions::StdExpansion::FwdTrans(), Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), Nektar::StdRegions::StdExpansion::GetTotPoints(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, and Vmath::Vcopy().

◆ v_FwdTrans()

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 220 of file PyrExp.cpp.

222  {
223  if(m_base[0]->Collocation() &&
224  m_base[1]->Collocation() &&
225  m_base[2]->Collocation())
226  {
227  Vmath::Vcopy(GetNcoeffs(),&inarray[0],1,&outarray[0],1);
228  }
229  else
230  {
231  v_IProductWRTBase(inarray,outarray);
232 
233  // get Mass matrix inverse
234  MatrixKey masskey(StdRegions::eInvMass,
235  DetShapeType(),*this);
236  DNekScalMatSharedPtr matsys = m_matrixManager[masskey];
237 
238  // copy inarray in case inarray == outarray
239  DNekVec in (m_ncoeffs,outarray);
240  DNekVec out(m_ncoeffs,outarray,eWrapper);
241 
242  out = (*matsys)*in;
243  }
244  }
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:278
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:124
NekVector< NekDouble > DNekVec
Definition: NekTypeDefs.hpp:48

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

◆ v_GenMatrix()

DNekMatSharedPtr Nektar::LocalRegions::PyrExp::v_GenMatrix ( const StdRegions::StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdPyrExp.

Definition at line 1073 of file PyrExp.cpp.

1074  {
1075  DNekMatSharedPtr returnval;
1076 
1077  switch(mkey.GetMatrixType())
1078  {
1085  returnval = Expansion3D::v_GenMatrix(mkey);
1086  break;
1087  default:
1088  returnval = StdPyrExp::v_GenMatrix(mkey);
1089  }
1090 
1091  return returnval;
1092  }
virtual DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey)

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

◆ v_GetCoord()

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 532 of file PyrExp.cpp.

534  {
535  int i;
536 
537  ASSERTL1(Lcoords[0] <= -1.0 && Lcoords[0] >= 1.0 &&
538  Lcoords[1] <= -1.0 && Lcoords[1] >= 1.0 &&
539  Lcoords[2] <= -1.0 && Lcoords[2] >= 1.0,
540  "Local coordinates are not in region [-1,1]");
541 
542  // m_geom->FillGeom(); // TODO: implement FillGeom()
543 
544  for(i = 0; i < m_geom->GetCoordim(); ++i)
545  {
546  coords[i] = m_geom->GetCoord(i,Lcoords);
547  }
548  }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:272

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

◆ v_GetCoordim()

int Nektar::LocalRegions::PyrExp::v_GetCoordim ( void  )
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 645 of file PyrExp.cpp.

646  {
647  return m_geom->GetCoordim();
648  }

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

◆ v_GetCoords()

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 550 of file PyrExp.cpp.

554  {
555  Expansion::v_GetCoords(coords_1, coords_2, coords_3);
556  }
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3)
Definition: Expansion.cpp:318

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

◆ v_GetLinStdExp()

StdRegions::StdExpansionSharedPtr Nektar::LocalRegions::PyrExp::v_GetLinStdExp ( void  ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 515 of file PyrExp.cpp.

516  {
517  LibUtilities::BasisKey bkey0(m_base[0]->GetBasisType(),
518  2, m_base[0]->GetPointsKey());
519  LibUtilities::BasisKey bkey1(m_base[1]->GetBasisType(),
520  2, m_base[1]->GetPointsKey());
521  LibUtilities::BasisKey bkey2(m_base[2]->GetBasisType(),
522  2, m_base[2]->GetPointsKey());
523 
525  ::AllocateSharedPtr( bkey0, bkey1, bkey2);
526  }

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

◆ v_GetLocMatrix()

DNekScalMatSharedPtr Nektar::LocalRegions::PyrExp::v_GetLocMatrix ( const MatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1105 of file PyrExp.cpp.

1106  {
1107  return m_matrixManager[mkey];
1108  }

References m_matrixManager.

◆ v_GetLocStaticCondMatrix()

DNekScalBlkMatSharedPtr Nektar::LocalRegions::PyrExp::v_GetLocStaticCondMatrix ( const MatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1110 of file PyrExp.cpp.

1111  {
1112  return m_staticCondMatrixManager[mkey];
1113  }

References m_staticCondMatrixManager.

◆ v_GetStdExp()

StdRegions::StdExpansionSharedPtr Nektar::LocalRegions::PyrExp::v_GetStdExp ( void  ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 507 of file PyrExp.cpp.

508  {
510  ::AllocateSharedPtr(m_base[0]->GetBasisKey(),
511  m_base[1]->GetBasisKey(),
512  m_base[2]->GetBasisKey());
513  }

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

◆ v_GetTracePhysMap()

void Nektar::LocalRegions::PyrExp::v_GetTracePhysMap ( const int  face,
Array< OneD, int > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 650 of file PyrExp.cpp.

652  {
653  int nquad0 = m_base[0]->GetNumPoints();
654  int nquad1 = m_base[1]->GetNumPoints();
655  int nquad2 = m_base[2]->GetNumPoints();
656 
657  int nq0 = 0;
658  int nq1 = 0;
659 
660  switch(face)
661  {
662  case 0:
663  nq0 = nquad0;
664  nq1 = nquad1;
665  if(outarray.size()!=nq0*nq1)
666  {
667  outarray = Array<OneD, int>(nq0*nq1);
668  }
669 
670  //Directions A and B positive
671  for(int i = 0; i < nquad0*nquad1; ++i)
672  {
673  outarray[i] = i;
674  }
675 
676  break;
677  case 1:
678  nq0 = nquad0;
679  nq1 = nquad2;
680  if(outarray.size()!=nq0*nq1)
681  {
682  outarray = Array<OneD, int>(nq0*nq1);
683  }
684 
685  //Direction A and B positive
686  for (int k=0; k<nquad2; k++)
687  {
688  for(int i = 0; i < nquad0; ++i)
689  {
690  outarray[k*nquad0+i] = (nquad0*nquad1*k)+i;
691  }
692  }
693 
694  break;
695  case 2:
696  nq0 = nquad1;
697  nq1 = nquad2;
698  if(outarray.size()!=nq0*nq1)
699  {
700  outarray = Array<OneD, int>(nq0*nq1);
701  }
702 
703  //Directions A and B positive
704  for(int j = 0; j < nquad1*nquad2; ++j)
705  {
706  outarray[j] = nquad0-1 + j*nquad0;
707 
708  }
709  break;
710  case 3:
711 
712  nq0 = nquad0;
713  nq1 = nquad2;
714  if(outarray.size()!=nq0*nq1)
715  {
716  outarray = Array<OneD, int>(nq0*nq1);
717  }
718 
719  //Direction A and B positive
720  for (int k=0; k<nquad2; k++)
721  {
722  for(int i = 0; i < nquad0; ++i)
723  {
724  outarray[k*nquad0+i] = nquad0*(nquad1-1) + (nquad0*nquad1*k)+i;
725  }
726  }
727  break;
728  case 4:
729  nq0 = nquad1;
730  nq1 = nquad2;
731 
732  if(outarray.size()!=nq0*nq1)
733  {
734  outarray = Array<OneD, int>(nq0*nq1);
735  }
736 
737  //Directions A and B positive
738  for(int j = 0; j < nquad1*nquad2; ++j)
739  {
740  outarray[j] = j*nquad0;
741 
742  }
743  break;
744  default:
745  ASSERTL0(false,"face value (> 4) is out of range");
746  break;
747  }
748  }

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

◆ v_Integral()

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 111 of file PyrExp.cpp.

112  {
113  int nquad0 = m_base[0]->GetNumPoints();
114  int nquad1 = m_base[1]->GetNumPoints();
115  int nquad2 = m_base[2]->GetNumPoints();
116  Array<OneD, const NekDouble> jac = m_metricinfo->GetJac(GetPointsKeys());
117  Array<OneD, NekDouble> tmp(nquad0*nquad1*nquad2);
118 
119  // multiply inarray with Jacobian
120  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
121  {
122  Vmath::Vmul(nquad0*nquad1*nquad2,&jac[0],1,(NekDouble*)&inarray[0],1, &tmp[0],1);
123  }
124  else
125  {
126  Vmath::Smul(nquad0*nquad1*nquad2,(NekDouble) jac[0], (NekDouble*)&inarray[0],1,&tmp[0],1);
127  }
128 
129  // call StdPyrExp version;
130  return StdPyrExp::v_Integral(tmp);
131  }

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

◆ v_IProductWRTBase()

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 278 of file PyrExp.cpp.

281  {
282  v_IProductWRTBase_SumFac(inarray, outarray);
283  }
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
Definition: PyrExp.cpp:285

References v_IProductWRTBase_SumFac().

Referenced by v_FwdTrans().

◆ v_IProductWRTBase_SumFac()

void Nektar::LocalRegions::PyrExp::v_IProductWRTBase_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
bool  multiplybyweights = true 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdPyrExp.

Definition at line 285 of file PyrExp.cpp.

289  {
290  const int nquad0 = m_base[0]->GetNumPoints();
291  const int nquad1 = m_base[1]->GetNumPoints();
292  const int nquad2 = m_base[2]->GetNumPoints();
293  const int order0 = m_base[0]->GetNumModes();
294  const int order1 = m_base[1]->GetNumModes();
295 
296  Array<OneD, NekDouble> wsp(order0*nquad2*(nquad1+order1));
297 
298  if(multiplybyweights)
299  {
300  Array<OneD, NekDouble> tmp(nquad0*nquad1*nquad2);
301 
302  MultiplyByQuadratureMetric(inarray, tmp);
303 
304  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
305  m_base[1]->GetBdata(),
306  m_base[2]->GetBdata(),
307  tmp,outarray,wsp,
308  true,true,true);
309  }
310  else
311  {
312  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
313  m_base[1]->GetBdata(),
314  m_base[2]->GetBdata(),
315  inarray,outarray,wsp,
316  true,true,true);
317  }
318  }
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)

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

Referenced by v_IProductWRTBase().

◆ v_IProductWRTDerivBase()

void Nektar::LocalRegions::PyrExp::v_IProductWRTDerivBase ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Calculates the inner product \( I_{pqr} = (u, \partial_{x_i} \phi_{pqr}) \).

The derivative of the basis functions is performed using the chain rule in order to incorporate the geometric factors. Assuming that the basis functions are a tensor product \(\phi_{pqr}(\eta_1,\eta_2,\eta_3) = \phi_1(\eta_1)\phi_2(\eta_2)\phi_3(\eta_3)\), this yields the result

\[ I_{pqr} = \sum_{j=1}^3 \left(u, \frac{\partial u}{\partial \eta_j} \frac{\partial \eta_j}{\partial x_i}\right) \]

In the pyramid element, we must also incorporate a second set of geometric factors which incorporate the collapsed co-ordinate system, so that

\[ \frac{\partial\eta_j}{\partial x_i} = \sum_{k=1}^3 \frac{\partial\eta_j}{\partial\xi_k}\frac{\partial\xi_k}{\partial x_i} \]

These derivatives can be found on p152 of Sherwin & Karniadakis.

Parameters
dirDirection in which to take the derivative.
inarrayThe function \( u \).
outarrayValue of the inner product.

Reimplemented from Nektar::StdRegions::StdPyrExp.

Definition at line 351 of file PyrExp.cpp.

355  {
356  v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
357  }
void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: PyrExp.cpp:359

References v_IProductWRTDerivBase_SumFac().

◆ v_IProductWRTDerivBase_SumFac()

void Nektar::LocalRegions::PyrExp::v_IProductWRTDerivBase_SumFac ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual
Parameters
inarrayFunction evaluated at physical collocation points.
outarrayInner product with respect to each basis function over the element.

Reimplemented from Nektar::StdRegions::StdPyrExp.

Definition at line 359 of file PyrExp.cpp.

363  {
364  const int nquad0 = m_base[0]->GetNumPoints();
365  const int nquad1 = m_base[1]->GetNumPoints();
366  const int nquad2 = m_base[2]->GetNumPoints();
367  const int order0 = m_base[0]->GetNumModes ();
368  const int order1 = m_base[1]->GetNumModes ();
369  const int nqtot = nquad0*nquad1*nquad2;
370 
371  Array<OneD, NekDouble> tmp1 (nqtot );
372  Array<OneD, NekDouble> tmp2 (nqtot );
373  Array<OneD, NekDouble> tmp3 (nqtot );
374  Array<OneD, NekDouble> tmp4 (nqtot );
375  Array<OneD, NekDouble> tmp6 (m_ncoeffs);
376  Array<OneD, NekDouble> wsp (std::max(nqtot,order0*nquad2*(nquad1+order1)));
377 
378  MultiplyByQuadratureMetric(inarray, tmp1);
379 
380  Array<OneD, Array<OneD, NekDouble>> tmp2D{3};
381  tmp2D[0] = tmp2;
382  tmp2D[1] = tmp3;
383  tmp2D[2] = tmp4;
384 
385  PyrExp::v_AlignVectorToCollapsedDir(dir, tmp1, tmp2D);
386 
387  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
388  m_base[1]->GetBdata (),
389  m_base[2]->GetBdata (),
390  tmp2,outarray,wsp,
391  false,true,true);
392 
393  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata (),
394  m_base[1]->GetDbdata(),
395  m_base[2]->GetBdata (),
396  tmp3,tmp6,wsp,
397  true,false,true);
398 
399  Vmath::Vadd(m_ncoeffs, tmp6, 1, outarray, 1, outarray, 1);
400 
401  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata (),
402  m_base[1]->GetBdata (),
403  m_base[2]->GetDbdata(),
404  tmp4,tmp6,wsp,
405  true,true,false);
406 
407  Vmath::Vadd(m_ncoeffs, tmp6, 1, outarray, 1, outarray, 1);
408  }
virtual void v_AlignVectorToCollapsedDir(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: PyrExp.cpp:410

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

Referenced by v_IProductWRTDerivBase().

◆ v_LaplacianMatrixOp_MatFree_Kernel()

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 1534 of file PyrExp.cpp.

1538  {
1539  // This implementation is only valid when there are no coefficients
1540  // associated to the Laplacian operator
1541  if (m_metrics.count(eMetricLaplacian00) == 0)
1542  {
1544  }
1545 
1546  int nquad0 = m_base[0]->GetNumPoints();
1547  int nquad1 = m_base[1]->GetNumPoints();
1548  int nq2 = m_base[2]->GetNumPoints();
1549  int nqtot = nquad0*nquad1*nq2;
1550 
1551  ASSERTL1(wsp.size() >= 6*nqtot,
1552  "Insufficient workspace size.");
1553  ASSERTL1(m_ncoeffs <= nqtot,
1554  "Workspace not set up for ncoeffs > nqtot");
1555 
1556  const Array<OneD, const NekDouble>& base0 = m_base[0]->GetBdata();
1557  const Array<OneD, const NekDouble>& base1 = m_base[1]->GetBdata();
1558  const Array<OneD, const NekDouble>& base2 = m_base[2]->GetBdata();
1559  const Array<OneD, const NekDouble>& dbase0 = m_base[0]->GetDbdata();
1560  const Array<OneD, const NekDouble>& dbase1 = m_base[1]->GetDbdata();
1561  const Array<OneD, const NekDouble>& dbase2 = m_base[2]->GetDbdata();
1562  const Array<OneD, const NekDouble>& metric00 = m_metrics[eMetricLaplacian00];
1563  const Array<OneD, const NekDouble>& metric01 = m_metrics[eMetricLaplacian01];
1564  const Array<OneD, const NekDouble>& metric02 = m_metrics[eMetricLaplacian02];
1565  const Array<OneD, const NekDouble>& metric11 = m_metrics[eMetricLaplacian11];
1566  const Array<OneD, const NekDouble>& metric12 = m_metrics[eMetricLaplacian12];
1567  const Array<OneD, const NekDouble>& metric22 = m_metrics[eMetricLaplacian22];
1568 
1569  // Allocate temporary storage
1570  Array<OneD,NekDouble> wsp0 (2*nqtot, wsp);
1571  Array<OneD,NekDouble> wsp1 ( nqtot, wsp+1*nqtot);
1572  Array<OneD,NekDouble> wsp2 ( nqtot, wsp+2*nqtot);
1573  Array<OneD,NekDouble> wsp3 ( nqtot, wsp+3*nqtot);
1574  Array<OneD,NekDouble> wsp4 ( nqtot, wsp+4*nqtot);
1575  Array<OneD,NekDouble> wsp5 ( nqtot, wsp+5*nqtot);
1576 
1577  // LAPLACIAN MATRIX OPERATION
1578  // wsp1 = du_dxi1 = D_xi1 * inarray = D_xi1 * u
1579  // wsp2 = du_dxi2 = D_xi2 * inarray = D_xi2 * u
1580  // wsp2 = du_dxi3 = D_xi3 * inarray = D_xi3 * u
1581  StdExpansion3D::PhysTensorDeriv(inarray,wsp0,wsp1,wsp2);
1582 
1583  // wsp0 = k = g0 * wsp1 + g1 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
1584  // wsp2 = l = g1 * wsp1 + g2 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
1585  // where g0, g1 and g2 are the metric terms set up in the GeomFactors class
1586  // especially for this purpose
1587  Vmath::Vvtvvtp(nqtot,&metric00[0],1,&wsp0[0],1,&metric01[0],1,&wsp1[0],1,&wsp3[0],1);
1588  Vmath::Vvtvp (nqtot,&metric02[0],1,&wsp2[0],1,&wsp3[0],1,&wsp3[0],1);
1589  Vmath::Vvtvvtp(nqtot,&metric01[0],1,&wsp0[0],1,&metric11[0],1,&wsp1[0],1,&wsp4[0],1);
1590  Vmath::Vvtvp (nqtot,&metric12[0],1,&wsp2[0],1,&wsp4[0],1,&wsp4[0],1);
1591  Vmath::Vvtvvtp(nqtot,&metric02[0],1,&wsp0[0],1,&metric12[0],1,&wsp1[0],1,&wsp5[0],1);
1592  Vmath::Vvtvp (nqtot,&metric22[0],1,&wsp2[0],1,&wsp5[0],1,&wsp5[0],1);
1593 
1594  // outarray = m = (D_xi1 * B)^T * k
1595  // wsp1 = n = (D_xi2 * B)^T * l
1596  IProductWRTBase_SumFacKernel(dbase0,base1,base2,wsp3,outarray,wsp0,false,true,true);
1597  IProductWRTBase_SumFacKernel(base0,dbase1,base2,wsp4,wsp2, wsp0,true,false,true);
1598  Vmath::Vadd(m_ncoeffs,wsp2.get(),1,outarray.get(),1,outarray.get(),1);
1599  IProductWRTBase_SumFacKernel(base0,base1,dbase2,wsp5,wsp2, wsp0,true,true,false);
1600  Vmath::Vadd(m_ncoeffs,wsp2.get(),1,outarray.get(),1,outarray.get(),1);
1601  }

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

◆ v_PhysDeriv()

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 138 of file PyrExp.cpp.

142  {
143  int nquad0 = m_base[0]->GetNumPoints();
144  int nquad1 = m_base[1]->GetNumPoints();
145  int nquad2 = m_base[2]->GetNumPoints();
146  Array<TwoD, const NekDouble> gmat =
147  m_metricinfo->GetDerivFactors(GetPointsKeys());
148  Array<OneD,NekDouble> diff0(nquad0*nquad1*nquad2);
149  Array<OneD,NekDouble> diff1(nquad0*nquad1*nquad2);
150  Array<OneD,NekDouble> diff2(nquad0*nquad1*nquad2);
151 
152  StdPyrExp::v_PhysDeriv(inarray, diff0, diff1, diff2);
153 
154  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
155  {
156  if(out_d0.size())
157  {
158  Vmath::Vmul (nquad0*nquad1*nquad2,&gmat[0][0],1,&diff0[0],1, &out_d0[0], 1);
159  Vmath::Vvtvp (nquad0*nquad1*nquad2,&gmat[1][0],1,&diff1[0],1, &out_d0[0], 1,&out_d0[0],1);
160  Vmath::Vvtvp (nquad0*nquad1*nquad2,&gmat[2][0],1,&diff2[0],1, &out_d0[0], 1,&out_d0[0],1);
161  }
162 
163  if(out_d1.size())
164  {
165  Vmath::Vmul (nquad0*nquad1*nquad2,&gmat[3][0],1,&diff0[0],1, &out_d1[0], 1);
166  Vmath::Vvtvp (nquad0*nquad1*nquad2,&gmat[4][0],1,&diff1[0],1, &out_d1[0], 1,&out_d1[0],1);
167  Vmath::Vvtvp (nquad0*nquad1*nquad2,&gmat[5][0],1,&diff2[0],1, &out_d1[0], 1,&out_d1[0],1);
168  }
169 
170  if(out_d2.size())
171  {
172  Vmath::Vmul (nquad0*nquad1*nquad2,&gmat[6][0],1,&diff0[0],1, &out_d2[0], 1);
173  Vmath::Vvtvp (nquad0*nquad1*nquad2,&gmat[7][0],1,&diff1[0],1, &out_d2[0], 1, &out_d2[0],1);
174  Vmath::Vvtvp (nquad0*nquad1*nquad2,&gmat[8][0],1,&diff2[0],1, &out_d2[0], 1, &out_d2[0],1);
175  }
176  }
177  else // regular geometry
178  {
179  if(out_d0.size())
180  {
181  Vmath::Smul (nquad0*nquad1*nquad2,gmat[0][0],&diff0[0],1, &out_d0[0], 1);
182  Blas::Daxpy (nquad0*nquad1*nquad2,gmat[1][0],&diff1[0],1, &out_d0[0], 1);
183  Blas::Daxpy (nquad0*nquad1*nquad2,gmat[2][0],&diff2[0],1, &out_d0[0], 1);
184  }
185 
186  if(out_d1.size())
187  {
188  Vmath::Smul (nquad0*nquad1*nquad2,gmat[3][0],&diff0[0],1, &out_d1[0], 1);
189  Blas::Daxpy (nquad0*nquad1*nquad2,gmat[4][0],&diff1[0],1, &out_d1[0], 1);
190  Blas::Daxpy (nquad0*nquad1*nquad2,gmat[5][0],&diff2[0],1, &out_d1[0], 1);
191  }
192 
193  if(out_d2.size())
194  {
195  Vmath::Smul (nquad0*nquad1*nquad2,gmat[6][0],&diff0[0],1, &out_d2[0], 1);
196  Blas::Daxpy (nquad0*nquad1*nquad2,gmat[7][0],&diff1[0],1, &out_d2[0], 1);
197  Blas::Daxpy (nquad0*nquad1*nquad2,gmat[8][0],&diff2[0],1, &out_d2[0], 1);
198  }
199  }
200  }
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
Definition: Blas.hpp:167

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

◆ v_PhysEvaluate()

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 627 of file PyrExp.cpp.

629  {
630  Array<OneD,NekDouble> Lcoord(3);
631 
632  ASSERTL0(m_geom,"m_geom not defined");
633 
634  //TODO: check GetLocCoords()
635  m_geom->GetLocCoords(coord, Lcoord);
636 
637  return StdPyrExp::v_PhysEvaluate(Lcoord, physvals);
638  }

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

◆ v_StdPhysEvaluate()

NekDouble Nektar::LocalRegions::PyrExp::v_StdPhysEvaluate ( const Array< OneD, const NekDouble > &  Lcoord,
const Array< OneD, const NekDouble > &  physvals 
)
protectedvirtual

Given the local cartesian coordinate Lcoord evaluate the value of physvals at this point by calling through to the StdExpansion method

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 619 of file PyrExp.cpp.

622  {
623  // Evaluate point in local coordinates.
624  return StdPyrExp::v_PhysEvaluate(Lcoord,physvals);
625  }

◆ v_SVVLaplacianFilter()

void Nektar::LocalRegions::PyrExp::v_SVVLaplacianFilter ( Array< OneD, NekDouble > &  array,
const StdRegions::StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdPyrExp.

Definition at line 1039 of file PyrExp.cpp.

1042  {
1043  int nq = GetTotPoints();
1044 
1045  // Calculate sqrt of the Jacobian
1046  Array<OneD, const NekDouble> jac =
1047  m_metricinfo->GetJac(GetPointsKeys());
1048  Array<OneD, NekDouble> sqrt_jac(nq);
1049  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1050  {
1051  Vmath::Vsqrt(nq,jac,1,sqrt_jac,1);
1052  }
1053  else
1054  {
1055  Vmath::Fill(nq,sqrt(jac[0]),sqrt_jac,1);
1056  }
1057 
1058  // Multiply array by sqrt(Jac)
1059  Vmath::Vmul(nq,sqrt_jac,1,array,1,array,1);
1060 
1061  // Apply std region filter
1062  StdPyrExp::v_SVVLaplacianFilter( array, mkey);
1063 
1064  // Divide by sqrt(Jac)
1065  Vmath::Vdiv(nq,array,1,sqrt_jac,1,array,1);
1066  }
void Vdiv(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:257

References Nektar::SpatialDomains::eDeformed, Vmath::Fill(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::GetTotPoints(), Nektar::LocalRegions::Expansion::m_metricinfo, tinysimd::sqrt(), Vmath::Vdiv(), Vmath::Vmul(), and Vmath::Vsqrt().

Member Data Documentation

◆ m_matrixManager

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

Definition at line 186 of file PyrExp.h.

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

◆ m_staticCondMatrixManager

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

Definition at line 187 of file PyrExp.h.

Referenced by v_DropLocStaticCondMatrix(), and v_GetLocStaticCondMatrix().