Nektar++
Public Member Functions | Protected Member Functions | Private Types | List of all members
Nektar::StdRegions::StdQuadExp Class Reference

#include <StdQuadExp.h>

Inheritance diagram for Nektar::StdRegions::StdQuadExp:
[legend]

Public Member Functions

 StdQuadExp ()
 
 StdQuadExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb)
 Constructor using BasisKey class for quadrature points and order definition. More...
 
 StdQuadExp (const StdQuadExp &T)
 Copy Constructor. More...
 
 ~StdQuadExp ()
 Destructor. More...
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion2D
 StdExpansion2D ()
 
 StdExpansion2D (int numcoeffs, const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb)
 
 StdExpansion2D (const StdExpansion2D &T)
 
virtual ~StdExpansion2D ()
 
void PhysTensorDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d0, Array< OneD, NekDouble > &outarray_d1)
 Calculate the 2D derivative in the local tensor/collapsed coordinate at the physical points. More...
 
NekDouble Integral (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &w0, const Array< OneD, const NekDouble > &w1)
 
void BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
 
void IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion
 StdExpansion ()
 Default Constructor. More...
 
 StdExpansion (const int numcoeffs, const int numbases, const LibUtilities::BasisKey &Ba=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bb=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bc=LibUtilities::NullBasisKey)
 Constructor. More...
 
 StdExpansion (const StdExpansion &T)
 Copy Constructor. More...
 
virtual ~StdExpansion ()
 Destructor. More...
 
int GetNumBases () const
 This function returns the number of 1D bases used in the expansion. More...
 
const Array< OneD, const LibUtilities::BasisSharedPtr > & GetBase () const
 This function gets the shared point to basis. More...
 
const LibUtilities::BasisSharedPtrGetBasis (int dir) const
 This function gets the shared point to basis in the dir direction. More...
 
int GetNcoeffs (void) const
 This function returns the total number of coefficients used in the expansion. More...
 
int GetTotPoints () const
 This function returns the total number of quadrature points used in the element. More...
 
LibUtilities::BasisType GetBasisType (const int dir) const
 This function returns the type of basis used in the dir direction. More...
 
int GetBasisNumModes (const int dir) const
 This function returns the number of expansion modes in the dir direction. More...
 
int EvalBasisNumModesMax (void) const
 This function returns the maximum number of expansion modes over all local directions. More...
 
LibUtilities::PointsType GetPointsType (const int dir) const
 This function returns the type of quadrature points used in the dir direction. More...
 
int GetNumPoints (const int dir) const
 This function returns the number of quadrature points in the dir direction. More...
 
const Array< OneD, const NekDouble > & GetPoints (const int dir) const
 This function returns a pointer to the array containing the quadrature points in dir direction. More...
 
int GetNverts () const
 This function returns the number of vertices of the expansion domain. More...
 
int GetNedges () const
 This function returns the number of edges of the expansion domain. More...
 
int GetEdgeNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th edge. More...
 
int GetTotalEdgeIntNcoeffs () const
 
int GetEdgeNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th edge. More...
 
int DetCartesianDirOfEdge (const int edge)
 
const LibUtilities::BasisKey DetEdgeBasisKey (const int i) const
 
const LibUtilities::BasisKey DetFaceBasisKey (const int i, const int k) const
 
int GetFaceNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th face. More...
 
int GetFaceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th face. More...
 
int GetFaceIntNcoeffs (const int i) const
 
int GetTotalFaceIntNcoeffs () const
 
int GetTraceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th edge/face. More...
 
LibUtilities::PointsKey GetFacePointsKey (const int i, const int j) const
 
int NumBndryCoeffs (void) const
 
int NumDGBndryCoeffs (void) const
 
LibUtilities::BasisType GetEdgeBasisType (const int i) const
 This function returns the type of expansion basis on the i-th edge. More...
 
const LibUtilities::PointsKey GetNodalPointsKey () const
 This function returns the type of expansion Nodal point type if defined. More...
 
int GetNfaces () const
 This function returns the number of faces of the expansion domain. More...
 
int GetNtrace () const
 Returns the number of trace elements connected to this element. More...
 
LibUtilities::ShapeType DetShapeType () const
 This function returns the shape of the expansion domain. More...
 
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)
 
IndexMapValuesSharedPtr GetIndexMap (const IndexMapKey &ikey)
 
const Array< OneD, const NekDouble > & GetPhysNormals (void)
 
void SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
virtual void SetUpPhysNormals (const int edge)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
void DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
StdRegions::Orientation GetForient (int face)
 
StdRegions::Orientation GetEorient (int edge)
 
void SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
 
void SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
int CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
NekDouble StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
int GetCoordim ()
 
void GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
void GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
int GetVertexMap (const int localVertexId, bool useCoeffPacking=false)
 
void GetEdgeInteriorMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
 
void GetFaceNumModes (const int fid, const Orientation faceOrient, int &numModes0, int &numModes1)
 
void GetFaceInteriorMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
 
void GetEdgeToElementMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int P=-1)
 
void GetFaceToElementMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int nummodesA=-1, int nummodesB=-1)
 
void GetEdgePhysVals (const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Extract the physical values along edge edge from inarray into outarray following the local edge orientation and point distribution defined by defined in EdgeExp. More...
 
void GetEdgePhysVals (const int edge, const std::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetTracePhysVals (const int edge, const std::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetVertexPhysVals (const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray)
 
void GetEdgeInterpVals (const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetEdgeQFactors (const int edge, Array< OneD, NekDouble > &outarray)
 Extract the metric factors to compute the contravariant fluxes along edge edge and stores them into outarray following the local edge orientation (i.e. anticlockwise convention). More...
 
void GetFacePhysVals (const int face, const std::shared_ptr< StdExpansion > &FaceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=eNoOrientation)
 
void GetEdgePhysMap (const int edge, Array< OneD, int > &outarray)
 
void GetFacePhysMap (const int face, Array< OneD, int > &outarray)
 
void MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
DNekMatSharedPtr CreateGeneralMatrix (const StdMatrixKey &mkey)
 this function generates the mass matrix \(\mathbf{M}[i][j] = \int \phi_i(\mathbf{x}) \phi_j(\mathbf{x}) d\mathbf{x}\) More...
 
void GeneralMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
void 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)
 
void AddRobinMassMatrix (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
void AddRobinEdgeContribution (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, Array< OneD, NekDouble > &coeffs)
 
NekDouble PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta. More...
 
virtual int v_GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
virtual const Array< OneD, const NekDouble > & v_GetPhysNormals (void)
 
virtual void v_SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
virtual void v_SetUpPhysNormals (const int edge)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual void v_DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual StdRegions::Orientation v_GetForient (int face)
 
virtual StdRegions::Orientation v_GetEorient (int edge)
 
NekDouble Linf (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( L_\infty\) error \( |\epsilon|_\infty = \max |u - u_{exact}|\) where \( u_{exact}\) is given by the array sol. More...
 
NekDouble L2 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( L_2\) error, \( | \epsilon |_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 dx \right]^{1/2} d\xi_1 \) where \( u_{exact}\) is given by the array sol. More...
 
NekDouble H1 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( H^1\) error, \( | \epsilon |^1_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 + \nabla(u - u_{exact})\cdot\nabla(u - u_{exact})\cdot dx \right]^{1/2} d\xi_1 \) where \( u_{exact}\) is given by the array sol. More...
 
const NormalVectorGetEdgeNormal (const int edge) const
 
void ComputeEdgeNormal (const int edge)
 
void NegateEdgeNormal (const int edge)
 
bool EdgeNormalNegated (const int edge)
 
void ComputeFaceNormal (const int face)
 
void NegateFaceNormal (const int face)
 
bool FaceNormalNegated (const int face)
 
void ComputeVertexNormal (const int vertex)
 
void NegateVertexNormal (const int vertex)
 
bool VertexNormalNegated (const int vertex)
 
const NormalVectorGetFaceNormal (const int face) const
 
const NormalVectorGetVertexNormal (const int vertex) const
 
const NormalVectorGetSurfaceNormal (const int id) const
 
const LibUtilities::PointsKeyVector GetPointsKeys () const
 
Array< OneD, unsigned int > GetEdgeInverseBoundaryMap (int eid)
 
Array< OneD, unsigned int > GetFaceInverseBoundaryMap (int fid, StdRegions::Orientation faceOrient=eNoOrientation, 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)
 
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)
 

Protected Member Functions

NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray)
 Integrates the specified function over the domain. 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=NullNekDouble1DArray)
 Calculate the derivative of the physical points. More...
 
virtual void v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Calculate the derivative of the physical points in a given direction. More...
 
virtual void v_StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
virtual void v_StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_FwdTrans_BndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1)
 
virtual void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Transform a given function from physical quadrature space to coefficient space. 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 and put into outarray. More...
 
virtual void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
virtual void v_IProductWRTBase_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1)
 
virtual void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTDerivBase_MatOp (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 
virtual void v_FillMode (const int mode, Array< OneD, NekDouble > &array)
 Fill outarray with mode mode of expansion. More...
 
virtual int v_GetNverts () const
 
virtual int v_GetNedges () const
 
virtual int v_GetEdgeNcoeffs (const int i) const
 
virtual int v_GetEdgeNumPoints (const int i) const
 
virtual int v_NumBndryCoeffs () const
 
virtual int v_NumDGBndryCoeffs () const
 
virtual LibUtilities::BasisType v_GetEdgeBasisType (const int i) const
 
virtual int v_DetCartesianDirOfEdge (const int edge)
 
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
virtual const LibUtilities::BasisKey v_DetEdgeBasisKey (const int i) const
 
virtual LibUtilities::ShapeType v_DetShapeType () const
 
virtual bool v_IsBoundaryInteriorExpansion ()
 
virtual void v_GetCoords (Array< OneD, NekDouble > &coords_0, Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2)
 
virtual void v_GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
virtual void v_GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
virtual int v_GetVertexMap (int localVertexId, bool useCoeffPacking=false)
 
virtual void v_GetEdgeInteriorMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
 
void v_GetEdgeToElementMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int P=-1)
 
virtual DNekMatSharedPtr v_GenMatrix (const StdMatrixKey &mkey)
 
virtual DNekMatSharedPtr v_CreateStdMatrix (const StdMatrixKey &mkey)
 
virtual void v_MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
virtual void v_ExponentialFilter (Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
 
virtual void v_ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void v_GeneralMatrixOp_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true)
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion2D
virtual NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
virtual NekDouble v_PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 
virtual void v_LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
 
virtual void v_HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
 
virtual int v_GetTraceNcoeffs (const int i) const
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion
DNekMatSharedPtr CreateStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr CreateStdStaticCondMatrix (const StdMatrixKey &mkey)
 Create the static condensation of a matrix when using a boundary interior decomposition. More...
 
IndexMapValuesSharedPtr CreateIndexMap (const IndexMapKey &ikey)
 Create an IndexMap which contains mapping information linking any specific element shape with either its boundaries, edges, faces, verteces, etc. More...
 
void BwdTrans_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void 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 (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
 
virtual NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 

Private Types

typedef std::shared_ptr< StdExpansion1DStdExpansion1DSharedPtr
 

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
 
LibUtilities::NekManager< IndexMapKey, IndexMapValues, IndexMapKey::opLessm_IndexMapManager
 

Detailed Description

Definition at line 50 of file StdQuadExp.h.

Member Typedef Documentation

◆ StdExpansion1DSharedPtr

Definition at line 53 of file StdQuadExp.h.

Constructor & Destructor Documentation

◆ StdQuadExp() [1/3]

Nektar::StdRegions::StdQuadExp::StdQuadExp ( )

Definition at line 50 of file StdQuadExp.cpp.

51  {
52  }

◆ StdQuadExp() [2/3]

Nektar::StdRegions::StdQuadExp::StdQuadExp ( const LibUtilities::BasisKey Ba,
const LibUtilities::BasisKey Bb 
)

Constructor using BasisKey class for quadrature points and order definition.

Definition at line 57 of file StdQuadExp.cpp.

58  :
59  StdExpansion (Ba.GetNumModes()*Bb.GetNumModes(),2,Ba,Bb),
60  StdExpansion2D(Ba.GetNumModes()*Bb.GetNumModes(),Ba,Bb)
61  {
62  }
StdExpansion()
Default Constructor.

◆ StdQuadExp() [3/3]

Nektar::StdRegions::StdQuadExp::StdQuadExp ( const StdQuadExp T)

Copy Constructor.

Definition at line 65 of file StdQuadExp.cpp.

65  :
66  StdExpansion(T),
68  {
69  }
StdExpansion()
Default Constructor.

◆ ~StdQuadExp()

Nektar::StdRegions::StdQuadExp::~StdQuadExp ( )

Destructor.

Definition at line 72 of file StdQuadExp.cpp.

73  {
74  }

Member Function Documentation

◆ v_BwdTrans()

void Nektar::StdRegions::StdQuadExp::v_BwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 156 of file StdQuadExp.cpp.

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

158  {
159  if(m_base[0]->Collocation() && m_base[1]->Collocation())
160  {
162  inarray, 1, outarray, 1);
163  }
164  else
165  {
166  StdQuadExp::v_BwdTrans_SumFac(inarray,outarray);
167  }
168  }
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:228
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdQuadExp.cpp:170
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064

◆ v_BwdTrans_SumFac()

void Nektar::StdRegions::StdQuadExp::v_BwdTrans_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 170 of file StdQuadExp.cpp.

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

Referenced by v_BwdTrans().

173  {
174  Array<OneD, NekDouble> wsp(m_base[0]->GetNumPoints()*
175  m_base[1]->GetNumModes());
176 
177  BwdTrans_SumFacKernel(m_base[0]->GetBdata(),
178  m_base[1]->GetBdata(),
179  inarray,outarray,wsp,true,true);
180  }
void BwdTrans_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:228
Array< OneD, LibUtilities::BasisSharedPtr > m_base

◆ v_BwdTrans_SumFacKernel()

void Nektar::StdRegions::StdQuadExp::v_BwdTrans_SumFacKernel ( const Array< OneD, const NekDouble > &  base0,
const Array< OneD, const NekDouble > &  base1,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
Array< OneD, NekDouble > &  wsp,
bool  doCheckCollDir0,
bool  doCheckCollDir1 
)
protectedvirtual

Implements Nektar::StdRegions::StdExpansion2D.

Definition at line 192 of file StdQuadExp.cpp.

References ASSERTL1, Blas::Dgemm(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, and Vmath::Vcopy().

200  {
201  int nquad0 = m_base[0]->GetNumPoints();
202  int nquad1 = m_base[1]->GetNumPoints();
203  int nmodes0 = m_base[0]->GetNumModes();
204  int nmodes1 = m_base[1]->GetNumModes();
205 
206  bool colldir0 = doCheckCollDir0?(m_base[0]->Collocation()):false;
207  bool colldir1 = doCheckCollDir1?(m_base[1]->Collocation()):false;
208 
209  if(colldir0 && colldir1)
210  {
211  Vmath::Vcopy(m_ncoeffs,inarray.get(),1,outarray.get(),1);
212  }
213  else if(colldir0)
214  {
215  Blas::Dgemm('N','T', nquad0, nquad1,nmodes1, 1.0, &inarray[0], nquad0,
216  base1.get(), nquad1, 0.0, &outarray[0], nquad0);
217  }
218  else if(colldir1)
219  {
220  Blas::Dgemm('N','N', nquad0,nmodes1,nmodes0,1.0, base0.get(),
221  nquad0, &inarray[0], nmodes0,0.0,&outarray[0], nquad0);
222  }
223  else
224  {
225  ASSERTL1(wsp.num_elements()>=nquad0*nmodes1,"Workspace size is not sufficient");
226 
227  // Those two calls correpsond to the operation
228  // out = B0*in*Transpose(B1);
229  Blas::Dgemm('N','N', nquad0,nmodes1,nmodes0,1.0, base0.get(),
230  nquad0, &inarray[0], nmodes0,0.0,&wsp[0], nquad0);
231  Blas::Dgemm('N','T', nquad0, nquad1,nmodes1, 1.0, &wsp[0], nquad0,
232  base1.get(), nquad1, 0.0, &outarray[0], nquad0);
233  }
234  }
static void Dgemm(const char &transa, const char &transb, const int &m, const int &n, const int &k, const double &alpha, const double *a, const int &lda, const double *b, const int &ldb, const double &beta, double *c, const int &ldc)
BLAS level 3: Matrix-matrix multiply C = A x B where A[m x n], B[n x k], C[m x k].
Definition: Blas.hpp:213
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064

◆ v_CalcNumberOfCoefficients()

int Nektar::StdRegions::StdQuadExp::v_CalcNumberOfCoefficients ( const std::vector< unsigned int > &  nummodes,
int &  modes_offset 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 737 of file StdQuadExp.cpp.

740  {
741  int nmodes = nummodes[modes_offset]*nummodes[modes_offset+1];
742  modes_offset += 2;
743 
744  return nmodes;
745  }

◆ v_CreateStdMatrix()

DNekMatSharedPtr Nektar::StdRegions::StdQuadExp::v_CreateStdMatrix ( const StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 1511 of file StdQuadExp.cpp.

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

1512  {
1513  return GenMatrix(mkey);
1514  }
DNekMatSharedPtr GenMatrix(const StdMatrixKey &mkey)

◆ v_DetCartesianDirOfEdge()

int Nektar::StdRegions::StdQuadExp::v_DetCartesianDirOfEdge ( const int  edge)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 674 of file StdQuadExp.cpp.

References ASSERTL2.

675  {
676  ASSERTL2((edge >= 0)&&(edge <= 3),"edge id is out of range");
677 
678  if((edge == 0)||(edge == 2))
679  {
680  return 0;
681  }
682  else
683  {
684  return 1;
685  }
686  }
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:274

◆ v_DetEdgeBasisKey()

const LibUtilities::BasisKey Nektar::StdRegions::StdQuadExp::v_DetEdgeBasisKey ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 688 of file StdQuadExp.cpp.

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

689  {
690  ASSERTL2((i >= 0)&&(i <= 3),"edge id is out of range");
691 
692  if((i == 0)||(i == 2))
693  {
694  return GetBasis(0)->GetBasisKey();
695  }
696  else
697  {
698  return GetBasis(1)->GetBasisKey();
699  }
700 
701  }
const LibUtilities::BasisSharedPtr & GetBasis(int dir) const
This function gets the shared point to basis in the dir direction.
Definition: StdExpansion.h:117
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:274

◆ v_DetShapeType()

LibUtilities::ShapeType Nektar::StdRegions::StdQuadExp::v_DetShapeType ( ) const
protectedvirtual

◆ v_ExponentialFilter()

void Nektar::StdRegions::StdQuadExp::v_ExponentialFilter ( Array< OneD, NekDouble > &  array,
const NekDouble  alpha,
const NekDouble  exponent,
const NekDouble  cutoff 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1644 of file StdQuadExp.cpp.

References Nektar::StdRegions::StdExpansion::BwdTrans(), Nektar::LibUtilities::eOrtho_A, Nektar::StdRegions::StdExpansion::FwdTrans(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::StdRegions::StdExpansion::m_base, and class_topology::P.

1649  {
1650  // Generate an orthogonal expansion
1651  int qa = m_base[0]->GetNumPoints();
1652  int qb = m_base[1]->GetNumPoints();
1653  int nmodesA = m_base[0]->GetNumModes();
1654  int nmodesB = m_base[1]->GetNumModes();
1655  int P = nmodesA - 1;
1656  int Q = nmodesB - 1;
1657 
1658  // Declare orthogonal basis.
1659  LibUtilities::PointsKey pa(qa,m_base[0]->GetPointsType());
1660  LibUtilities::PointsKey pb(qb,m_base[1]->GetPointsType());
1661 
1662  LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A, nmodesA, pa);
1663  LibUtilities::BasisKey Bb(LibUtilities::eOrtho_A, nmodesB, pb);
1664  StdQuadExp OrthoExp(Ba,Bb);
1665 
1666  // Cutoff
1667  int Pcut = cutoff*P;
1668  int Qcut = cutoff*Q;
1669 
1670  // Project onto orthogonal space.
1671  Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
1672  OrthoExp.FwdTrans(array,orthocoeffs);
1673 
1674  //
1675  NekDouble fac, fac1, fac2;
1676  for(int i = 0; i < nmodesA; ++i)
1677  {
1678  for(int j = 0; j < nmodesB; ++j)
1679  {
1680  //to filter out only the "high-modes"
1681  if(i > Pcut || j > Qcut)
1682  {
1683  fac1 = (NekDouble) (i - Pcut)/( (NekDouble)(P - Pcut) );
1684  fac2 = (NekDouble) (j - Qcut)/( (NekDouble)(Q - Qcut) );
1685  fac = max(fac1, fac2);
1686  fac = pow(fac, exponent);
1687  orthocoeffs[i*nmodesB+j] *= exp(-alpha*fac);
1688  }
1689  }
1690  }
1691 
1692  // backward transform to physical space
1693  OrthoExp.BwdTrans(orthocoeffs,array);
1694  }
Principle Orthogonal Functions .
Definition: BasisType.h:45
double NekDouble
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:215
Array< OneD, LibUtilities::BasisSharedPtr > m_base

◆ v_FillMode()

void Nektar::StdRegions::StdQuadExp::v_FillMode ( const int  mode,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Fill outarray with mode mode of expansion.

Note for quadrilateral expansions _base0 modes run fastest

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 586 of file StdQuadExp.cpp.

References ASSERTL2, Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Vmath::Vcopy(), and Vmath::Vmul().

588  {
589  int i;
590  int nquad0 = m_base[0]->GetNumPoints();
591  int nquad1 = m_base[1]->GetNumPoints();
592  Array<OneD, const NekDouble> base0 = m_base[0]->GetBdata();
593  Array<OneD, const NekDouble> base1 = m_base[1]->GetBdata();
594  int btmp0 = m_base[0]->GetNumModes();
595  int mode0 = mode%btmp0;
596  int mode1 = mode/btmp0;
597 
598 
599  ASSERTL2(mode1 == (int)floor((1.0*mode)/btmp0),
600  "Integer Truncation not Equiv to Floor");
601 
602  ASSERTL2(m_ncoeffs <= mode,
603  "calling argument mode is larger than total expansion order");
604 
605  for(i = 0; i < nquad1; ++i)
606  {
607  Vmath::Vcopy(nquad0,(NekDouble *)(base0.get() + mode0*nquad0),
608  1, &outarray[0]+i*nquad0,1);
609  }
610 
611  for(i = 0; i < nquad0; ++i)
612  {
613  Vmath::Vmul(nquad1,(NekDouble *)(base1.get() + mode1*nquad1),1,
614  &outarray[0]+i,nquad0,&outarray[0]+i,nquad0);
615  }
616  }
double NekDouble
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:274
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
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:186

◆ v_FwdTrans()

void Nektar::StdRegions::StdQuadExp::v_FwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Transform a given function from physical quadrature space to coefficient space.

See also
StdExpansion::FwdTrans

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 237 of file StdQuadExp.cpp.

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

239  {
240  if((m_base[0]->Collocation())&&(m_base[1]->Collocation()))
241  {
242  Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
243  }
244  else
245  {
246  StdQuadExp::v_IProductWRTBase(inarray,outarray);
247 
248  // get Mass matrix inverse
249  StdMatrixKey masskey(eInvMass,DetShapeType(),*this);
250  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
251 
252  // copy inarray in case inarray == outarray
253  NekVector<NekDouble> in(m_ncoeffs,outarray,eCopy);
254  NekVector<NekDouble> out(m_ncoeffs,outarray,eWrapper);
255 
256  out = (*matsys)*in;
257  }
258  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:469
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:714
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 and put into outarray...
Definition: StdQuadExp.cpp:385
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064

◆ v_FwdTrans_BndConstrained()

void Nektar::StdRegions::StdQuadExp::v_FwdTrans_BndConstrained ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 260 of file StdQuadExp.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::StdRegions::StdExpansion::DetShapeType(), Blas::Dgemv(), Nektar::StdRegions::eForwards, Nektar::StdRegions::eMass, Nektar::StdRegions::StdExpansion::GetEdgeToElementMap(), Nektar::StdRegions::StdExpansion::GetInteriorMap(), Nektar::StdRegions::StdExpansion::IProductWRTBase(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::StdRegions::StdExpansion::m_stdStaticCondMatrixManager, Nektar::StdRegions::StdExpansion::MassMatrixOp(), Nektar::StdRegions::StdExpansion::NumBndryCoeffs(), Nektar::rhs, sign, Vmath::Vcopy(), and Vmath::Vsub().

263  {
264  if((m_base[0]->Collocation())&&(m_base[1]->Collocation()))
265  {
266  Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
267  }
268  else
269  {
270  int i,j;
271  int npoints[2] = {m_base[0]->GetNumPoints(),
272  m_base[1]->GetNumPoints()};
273  int nmodes[2] = {m_base[0]->GetNumModes(),
274  m_base[1]->GetNumModes()};
275 
276  fill(outarray.get(), outarray.get()+m_ncoeffs, 0.0 );
277 
278  Array<OneD, NekDouble> physEdge[4];
279  Array<OneD, NekDouble> coeffEdge[4];
280  for(i = 0; i < 4; i++)
281  {
282  physEdge[i] = Array<OneD, NekDouble>(npoints[i%2]);
283  coeffEdge[i] = Array<OneD, NekDouble>(nmodes[i%2]);
284  }
285 
286  for(i = 0; i < npoints[0]; i++)
287  {
288  physEdge[0][i] = inarray[i];
289  physEdge[2][i] = inarray[npoints[0]*npoints[1]-1-i];
290  }
291 
292  for(i = 0; i < npoints[1]; i++)
293  {
294  physEdge[1][i] = inarray[npoints[0]-1+i*npoints[0]];
295  physEdge[3][i] = inarray[(npoints[1]-1)*npoints[0]-i*npoints[0]];
296  }
297 
300 
301  Array<OneD, unsigned int> mapArray;
302  Array<OneD, int> signArray;
303  NekDouble sign;
304 
305  for(i = 0; i < 4; i++)
306  {
307  segexp[i%2]->FwdTrans_BndConstrained(physEdge[i],coeffEdge[i]);
308 
309  GetEdgeToElementMap(i,eForwards,mapArray,signArray);
310  for(j=0; j < nmodes[i%2]; j++)
311  {
312  sign = (NekDouble) signArray[j];
313  outarray[ mapArray[j] ] = sign * coeffEdge[i][j];
314  }
315  }
316 
317  Array<OneD, NekDouble> tmp0(m_ncoeffs);
318  Array<OneD, NekDouble> tmp1(m_ncoeffs);
319 
320  StdMatrixKey masskey(eMass,DetShapeType(),*this);
321  MassMatrixOp(outarray,tmp0,masskey);
322  IProductWRTBase(inarray,tmp1);
323 
324  Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
325 
326  // get Mass matrix inverse (only of interior DOF)
327  // use block (1,1) of the static condensed system
328  // note: this block alreay contains the inverse matrix
329  DNekMatSharedPtr matsys = (m_stdStaticCondMatrixManager[masskey])->GetBlock(1,1);
330 
331  int nBoundaryDofs = NumBndryCoeffs();
332  int nInteriorDofs = m_ncoeffs - nBoundaryDofs;
333 
334 
335  Array<OneD, NekDouble> rhs(nInteriorDofs);
336  Array<OneD, NekDouble> result(nInteriorDofs);
337 
338  GetInteriorMap(mapArray);
339 
340  for(i = 0; i < nInteriorDofs; i++)
341  {
342  rhs[i] = tmp1[ mapArray[i] ];
343  }
344 
345  Blas::Dgemv('N',nInteriorDofs,nInteriorDofs,1.0, &(matsys->GetPtr())[0],
346  nInteriorDofs,rhs.get(),1,0.0,result.get(),1);
347 
348  for(i = 0; i < nInteriorDofs; i++)
349  {
350  outarray[ mapArray[i] ] = result[i];
351  }
352  }
353 
354  }
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:974
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:469
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:16
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 expa...
Definition: StdExpansion.h:634
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
LibUtilities::NekManager< StdMatrixKey, DNekBlkMat, StdMatrixKey::opLess > m_stdStaticCondMatrixManager
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:817
double NekDouble
static void Dgemv(const char &trans, const int &m, const int &n, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = A x where A[m x n].
Definition: Blas.hpp:168
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:346
void GetEdgeToElementMap(const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int P=-1)
Definition: StdExpansion.h:849
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs
std::shared_ptr< StdSegExp > StdSegExpSharedPtr
Definition: StdSegExp.h:46
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064

◆ v_GeneralMatrixOp_MatOp()

void Nektar::StdRegions::StdQuadExp::v_GeneralMatrixOp_MatOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdMatrixKey mkey 
)
protected

Definition at line 1520 of file StdQuadExp.cpp.

References Blas::Dgemv(), Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::StdRegions::StdExpansion::m_stdMatrixManager, and Vmath::Vcopy().

1524  {
1525 
1527 
1528  if(inarray.get() == outarray.get())
1529  {
1530  Array<OneD,NekDouble> tmp(m_ncoeffs);
1531  Vmath::Vcopy(m_ncoeffs,inarray.get(),1,tmp.get(),1);
1532 
1533  Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, 1.0, mat->GetPtr().get(),
1534  m_ncoeffs, tmp.get(), 1, 0.0, outarray.get(), 1);
1535  }
1536  else
1537  {
1538  Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, 1.0, mat->GetPtr().get(),
1539  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
1540  }
1541  }
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
static void Dgemv(const char &trans, const int &m, const int &n, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = A x where A[m x n].
Definition: Blas.hpp:168
LibUtilities::NekManager< StdMatrixKey, DNekMat, StdMatrixKey::opLess > m_stdMatrixManager
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064

◆ v_GenMatrix()

DNekMatSharedPtr Nektar::StdRegions::StdQuadExp::v_GenMatrix ( const StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 1375 of file StdQuadExp.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::LibUtilities::BasisManager(), Nektar::StdRegions::StdMatrixKey::ConstFactorExists(), Nektar::StdRegions::StdExpansion::CreateGeneralMatrix(), Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::StdRegions::eFactorConst, Nektar::StdRegions::eFactorGaussEdge, Nektar::LibUtilities::eFourier, Nektar::StdRegions::eFwdTrans, Nektar::LibUtilities::eGauss_Lagrange, Nektar::StdRegions::eGaussDG, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::StdRegions::eInvMass, Nektar::StdRegions::eIProductWRTBase, Nektar::StdRegions::eMass, Nektar::StdRegions::ePhysInterpToEquiSpaced, Nektar::StdRegions::StdExpansion::GetBasisNumModes(), Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdMatrixKey::GetConstFactors(), Nektar::StdRegions::StdMatrixKey::GetMatrixType(), Nektar::LibUtilities::StdQuadData::getNumberOfCoefficients(), Nektar::StdRegions::StdExpansion::GetStdMatrix(), Nektar::StdRegions::StdExpansion::LocCoordToLocCollapsed(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Vmath::Smul(), and Vmath::Vcopy().

1376  {
1377  int i;
1378  int order0 = GetBasisNumModes(0);
1379  int order1 = GetBasisNumModes(1);
1380  MatrixType mtype = mkey.GetMatrixType();
1381 
1382  DNekMatSharedPtr Mat;
1383 
1384  switch(mtype)
1385  {
1387  {
1388  int nq0 = m_base[0]->GetNumPoints();
1389  int nq1 = m_base[1]->GetNumPoints();
1390  int nq;
1391 
1392  // take definition from key
1393  if(mkey.ConstFactorExists(eFactorConst))
1394  {
1395  nq = (int) mkey.GetConstFactor(eFactorConst);
1396  }
1397  else
1398  {
1399  nq = max(nq0,nq1);
1400  }
1401 
1402  int neq = LibUtilities::StdQuadData::
1403  getNumberOfCoefficients(nq, nq);
1404  Array<OneD, Array<OneD, NekDouble> > coords(neq);
1405  Array<OneD, NekDouble> coll (2);
1406  Array<OneD, DNekMatSharedPtr> I (2);
1407  Array<OneD, NekDouble> tmp (nq0);
1408 
1410  AllocateSharedPtr(neq, nq0 * nq1);
1411  int cnt = 0;
1412 
1413  for(int i = 0; i < nq; ++i)
1414  {
1415  for(int j = 0; j < nq; ++j,++cnt)
1416  {
1417  coords[cnt] = Array<OneD, NekDouble>(2);
1418  coords[cnt][0] = -1.0 + 2*j/(NekDouble)(nq-1);
1419  coords[cnt][1] = -1.0 + 2*i/(NekDouble)(nq-1);
1420  }
1421  }
1422 
1423  for(int i = 0; i < neq; ++i)
1424  {
1425  LocCoordToLocCollapsed(coords[i],coll);
1426 
1427  I[0] = m_base[0]->GetI(coll);
1428  I[1] = m_base[1]->GetI(coll+1);
1429 
1430  // interpolate first coordinate direction
1431  for (int j = 0; j < nq1; ++j)
1432  {
1433  NekDouble fac = (I[1]->GetPtr())[j];
1434  Vmath::Smul(nq0,fac,I[0]->GetPtr(),1,tmp,1);
1435 
1436  Vmath::Vcopy(nq0, &tmp[0], 1,
1437  Mat->GetRawPtr()+j*nq0*neq+i,neq);
1438  }
1439 
1440  }
1441  break;
1442  }
1443  case eMass:
1444  {
1446  // For Fourier basis set the imaginary component of mean mode
1447  // to have a unit diagonal component in mass matrix
1449  {
1450  for(i = 0; i < order1; ++i)
1451  {
1452  (*Mat)(order0*i+1,i*order0+1) = 1.0;
1453  }
1454  }
1455 
1457  {
1458  for(i = 0; i < order0; ++i)
1459  {
1460  (*Mat)(order0+i ,order0+i) = 1.0;
1461  }
1462  }
1463  break;
1464  }
1465  case eFwdTrans:
1466  {
1468  StdMatrixKey iprodkey(eIProductWRTBase,DetShapeType(),*this);
1469  DNekMat &Iprod = *GetStdMatrix(iprodkey);
1470  StdMatrixKey imasskey(eInvMass,DetShapeType(),*this);
1471  DNekMat &Imass = *GetStdMatrix(imasskey);
1472 
1473  (*Mat) = Imass*Iprod;
1474  break;
1475  }
1476  case eGaussDG:
1477  {
1478  ConstFactorMap factors = mkey.GetConstFactors();
1479 
1480  int edge = (int)factors[StdRegions::eFactorGaussEdge];
1481  int dir = (edge + 1) % 2;
1482  int nCoeffs = m_base[dir]->GetNumModes();
1483 
1484  const LibUtilities::PointsKey BS_p(
1486  const LibUtilities::BasisKey BS_k(
1487  LibUtilities::eGauss_Lagrange, nCoeffs, BS_p);
1488 
1489  Array<OneD, NekDouble> coords(1, 0.0);
1490  coords[0] = (edge == 0 || edge == 3) ? -1.0 : 1.0;
1491 
1494  DNekMatSharedPtr m_Ix = basis->GetI(coords);
1495 
1497  1.0, nCoeffs);
1498  Vmath::Vcopy(nCoeffs, m_Ix->GetPtr(), 1, Mat->GetPtr(), 1);
1499  break;
1500  }
1501  default:
1502  {
1504  break;
1505  }
1506  }
1507 
1508  return Mat;
1509  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:469
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:177
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:55
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:294
Fourier Expansion .
Definition: BasisType.h:53
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
std::shared_ptr< Basis > BasisSharedPtr
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:714
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:48
BasisManagerT & BasisManager(void)
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:216
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
NekMatrix< NekDouble, StandardMatrixTag > DNekMat
Definition: NekTypeDefs.hpp:51
double NekDouble
int getNumberOfCoefficients(int Na, int Nb)
Definition: ShapeType.hpp:137
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:164
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064

◆ v_GetBoundaryMap()

void Nektar::StdRegions::StdQuadExp::v_GetBoundaryMap ( Array< OneD, unsigned int > &  outarray)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 786 of file StdQuadExp.cpp.

References ASSERTL0, Nektar::LibUtilities::eGauss_Lagrange, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), Nektar::StdRegions::StdExpansion::m_base, and Nektar::StdRegions::StdExpansion::NumBndryCoeffs().

787  {
788  int i;
789  int cnt=0;
790  int nummodes0, nummodes1;
791  int value1 = 0, value2 = 0;
792  if(outarray.num_elements()!=NumBndryCoeffs())
793  {
794  outarray = Array<OneD, unsigned int>(NumBndryCoeffs());
795  }
796 
797  nummodes0 = m_base[0]->GetNumModes();
798  nummodes1 = m_base[1]->GetNumModes();
799 
800  const LibUtilities::BasisType Btype0 = GetBasisType(0);
801  const LibUtilities::BasisType Btype1 = GetBasisType(1);
802 
803  switch(Btype1)
804  {
807  value1 = nummodes0;
808  break;
810  value1 = 2*nummodes0;
811  break;
812  default:
813  ASSERTL0(0,"Mapping array is not defined for this expansion");
814  break;
815  }
816 
817  for(i = 0; i < value1; i++)
818  {
819  outarray[i]=i;
820  }
821  cnt=value1;
822 
823  switch(Btype0)
824  {
827  value2 = value1+nummodes0-1;
828  break;
830  value2 = value1+1;
831  break;
832  default:
833  ASSERTL0(0,"Mapping array is not defined for this expansion");
834  break;
835  }
836 
837  for(i = 0; i < nummodes1-2; i++)
838  {
839  outarray[cnt++]=value1+i*nummodes0;
840  outarray[cnt++]=value2+i*nummodes0;
841  }
842 
843 
845  {
846  for(i = nummodes0*(nummodes1-1);i < GetNcoeffs(); i++)
847  {
848  outarray[cnt++] = i;
849  }
850  }
851  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
Principle Modified Functions .
Definition: BasisType.h:48
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:55
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:130
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:164
Lagrange for SEM basis .
Definition: BasisType.h:54
Array< OneD, LibUtilities::BasisSharedPtr > m_base

◆ v_GetCoords()

void Nektar::StdRegions::StdQuadExp::v_GetCoords ( Array< OneD, NekDouble > &  coords_0,
Array< OneD, NekDouble > &  coords_1,
Array< OneD, NekDouble > &  coords_2 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 764 of file StdQuadExp.cpp.

References Blas::Dcopy(), Vmath::Fill(), Nektar::StdRegions::StdExpansion::GetNumPoints(), and Nektar::StdRegions::StdExpansion::m_base.

767  {
768  boost::ignore_unused(coords_2);
769  Array<OneD, const NekDouble> z0 = m_base[0]->GetZ();
770  Array<OneD, const NekDouble> z1 = m_base[1]->GetZ();
771  int nq0 = GetNumPoints(0);
772  int nq1 = GetNumPoints(1);
773  int i;
774 
775  for(i = 0; i < nq1; ++i)
776  {
777  Blas::Dcopy(nq0,z0.get(), 1,&coords_0[0] + i*nq0,1);
778  Vmath::Fill(nq0,z1[i],&coords_1[0] + i*nq0,1);
779  }
780  }
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:228
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45
Array< OneD, LibUtilities::BasisSharedPtr > m_base
static void Dcopy(const int &n, const double *x, const int &incx, double *y, const int &incy)
BLAS level 1: Copy x to y.
Definition: Blas.hpp:103

◆ v_GetEdgeBasisType()

LibUtilities::BasisType Nektar::StdRegions::StdQuadExp::v_GetEdgeBasisType ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 660 of file StdQuadExp.cpp.

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

661  {
662  ASSERTL2((i >= 0)&&(i <= 3),"edge id is out of range");
663 
664  if((i == 0)||(i == 2))
665  {
666  return GetBasisType(0);
667  }
668  else
669  {
670  return GetBasisType(1);
671  }
672  }
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:274
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:164

◆ v_GetEdgeInteriorMap()

void Nektar::StdRegions::StdQuadExp::v_GetEdgeInteriorMap ( const int  eid,
const Orientation  edgeOrient,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1014 of file StdQuadExp.cpp.

References ASSERTL0, Nektar::StdRegions::eBackwards, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetEdgeBasisType(), Nektar::StdRegions::StdExpansion::GetEdgeNcoeffs(), and Nektar::StdRegions::StdExpansion::m_base.

1018  {
1019  int i;
1020  const int nummodes0 = m_base[0]->GetNumModes();
1021  const int nummodes1 = m_base[1]->GetNumModes();
1022  const int nEdgeIntCoeffs = GetEdgeNcoeffs(eid)-2;
1023  const LibUtilities::BasisType bType = GetEdgeBasisType(eid);
1024 
1025  if(maparray.num_elements() != nEdgeIntCoeffs)
1026  {
1027  maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
1028  }
1029 
1030  if(signarray.num_elements() != nEdgeIntCoeffs)
1031  {
1032  signarray = Array<OneD, int>(nEdgeIntCoeffs,1);
1033  }
1034  else
1035  {
1036  fill( signarray.get() , signarray.get()+nEdgeIntCoeffs, 1 );
1037  }
1038 
1039  if(bType == LibUtilities::eModified_A)
1040  {
1041  switch(eid)
1042  {
1043  case 0:
1044  {
1045  for(i = 0; i < nEdgeIntCoeffs; i++)
1046  {
1047  maparray[i] = i+2;
1048  }
1049 
1050  if(edgeOrient==eBackwards)
1051  {
1052  for(i = 1; i < nEdgeIntCoeffs; i+=2)
1053  {
1054  signarray[i] = -1;
1055  }
1056  }
1057  }
1058  break;
1059  case 1:
1060  {
1061  for(i = 0; i < nEdgeIntCoeffs; i++)
1062  {
1063  maparray[i] = (i+2)*nummodes0 + 1;
1064  }
1065 
1066  if(edgeOrient==eBackwards)
1067  {
1068  for(i = 1; i < nEdgeIntCoeffs; i+=2)
1069  {
1070  signarray[i] = -1;
1071  }
1072  }
1073  }
1074  break;
1075  case 2:
1076  {
1077  for(i = 0; i < nEdgeIntCoeffs; i++)
1078  {
1079  maparray[i] = nummodes0+i+2;
1080  }
1081 
1082  if(edgeOrient==eBackwards)
1083  {
1084  for(i = 1; i < nEdgeIntCoeffs; i+=2)
1085  {
1086  signarray[i] = -1;
1087  }
1088  }
1089  }
1090  break;
1091  case 3:
1092  {
1093  for(i = 0; i < nEdgeIntCoeffs; i++)
1094  {
1095  maparray[i] = (i+2)*nummodes0;
1096  }
1097 
1098  if(edgeOrient==eBackwards)
1099  {
1100  for(i = 1; i < nEdgeIntCoeffs; i+=2)
1101  {
1102  signarray[i] = -1;
1103  }
1104  }
1105  }
1106  break;
1107  default:
1108  ASSERTL0(false,"eid must be between 0 and 3");
1109  break;
1110  }
1111  }
1112  else if(bType == LibUtilities::eGLL_Lagrange)
1113  {
1114  switch(eid)
1115  {
1116  case 0:
1117  {
1118  for(i = 0; i < nEdgeIntCoeffs; i++)
1119  {
1120  maparray[i] = i+1;
1121  }
1122  }
1123  break;
1124  case 1:
1125  {
1126  for(i = 0; i < nEdgeIntCoeffs; i++)
1127  {
1128  maparray[i] = (i+2)*nummodes0 - 1;
1129  }
1130  }
1131  break;
1132  case 2:
1133  {
1134  for(i = 0; i < nEdgeIntCoeffs; i++)
1135  {
1136  maparray[i] = nummodes0*(nummodes1-1) + i + 1;
1137  }
1138  }
1139  break;
1140  case 3:
1141  {
1142  for(i = 0; i < nEdgeIntCoeffs; i++)
1143  {
1144  maparray[i] = nummodes0 * (i+1);
1145  }
1146  }
1147  break;
1148  default:
1149  ASSERTL0(false,"eid must be between 0 and 3");
1150  break;
1151  }
1152  if(edgeOrient == eBackwards)
1153  {
1154  reverse( maparray.get() , maparray.get()+nEdgeIntCoeffs );
1155  }
1156  }
1157  else
1158  {
1159  ASSERTL0(false,"Mapping not defined for this type of basis");
1160  }
1161 
1162  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
Principle Modified Functions .
Definition: BasisType.h:48
int GetEdgeNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th edge.
Definition: StdExpansion.h:286
Lagrange for SEM basis .
Definition: BasisType.h:54
Array< OneD, LibUtilities::BasisSharedPtr > m_base
LibUtilities::BasisType GetEdgeBasisType(const int i) const
This function returns the type of expansion basis on the i-th edge.
Definition: StdExpansion.h:412

◆ v_GetEdgeNcoeffs()

int Nektar::StdRegions::StdQuadExp::v_GetEdgeNcoeffs ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 632 of file StdQuadExp.cpp.

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

633  {
634  ASSERTL2((i >= 0)&&(i <= 3),"edge id is out of range");
635 
636  if((i == 0)||(i == 2))
637  {
638  return GetBasisNumModes(0);
639  }
640  else
641  {
642  return GetBasisNumModes(1);
643  }
644  }
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:177
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:274

◆ v_GetEdgeNumPoints()

int Nektar::StdRegions::StdQuadExp::v_GetEdgeNumPoints ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 646 of file StdQuadExp.cpp.

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

647  {
648  ASSERTL2((i >= 0)&&(i <= 3),"edge id is out of range");
649 
650  if((i == 0)||(i == 2))
651  {
652  return GetNumPoints(0);
653  }
654  else
655  {
656  return GetNumPoints(1);
657  }
658  }
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:228
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:274

◆ v_GetEdgeToElementMap()

void Nektar::StdRegions::StdQuadExp::v_GetEdgeToElementMap ( const int  eid,
const Orientation  edgeOrient,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray,
int  P = -1 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1164 of file StdQuadExp.cpp.

References ASSERTL0, Nektar::StdRegions::eBackwards, Nektar::LibUtilities::eGauss_Lagrange, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetEdgeBasisType(), Nektar::StdRegions::StdExpansion::m_base, and class_topology::P.

1170  {
1171  int i;
1172  int numModes=0;
1173  int order0 = m_base[0]->GetNumModes();
1174  int order1 = m_base[1]->GetNumModes();
1175 
1176  switch (eid)
1177  {
1178  case 0:
1179  case 2:
1180  numModes = order0;
1181  break;
1182  case 1:
1183  case 3:
1184  numModes = order1;
1185  break;
1186  default:
1187  ASSERTL0(false,"eid must be between 0 and 3");
1188  }
1189 
1190  bool checkForZeroedModes = false;
1191  if (P == -1)
1192  {
1193  P = numModes;
1194  }
1195  else if(P != numModes)
1196  {
1197  checkForZeroedModes = true;
1198  }
1199  const LibUtilities::BasisType bType = GetEdgeBasisType(eid);
1200 
1201 
1202  if (maparray.num_elements() != P)
1203  {
1204  maparray = Array<OneD, unsigned int>(P);
1205  }
1206 
1207  if(signarray.num_elements() != P)
1208  {
1209  signarray = Array<OneD, int>(P, 1);
1210  }
1211  else
1212  {
1213  fill(signarray.get(), signarray.get()+P, 1);
1214  }
1215 
1216  if (bType == LibUtilities::eModified_A)
1217  {
1218  switch (eid)
1219  {
1220  case 0:
1221  {
1222  for (i = 0; i < P; i++)
1223  {
1224  maparray[i] = i;
1225  }
1226 
1227  if (edgeOrient == eBackwards)
1228  {
1229  swap(maparray[0], maparray[1]);
1230 
1231  for(i = 3; i < P; i+=2)
1232  {
1233  signarray[i] = -1;
1234  }
1235  }
1236  }
1237  break;
1238  case 1:
1239  {
1240  for (i = 0; i < P; i++)
1241  {
1242  maparray[i] = i*order0 + 1;
1243  }
1244 
1245  if (edgeOrient == eBackwards)
1246  {
1247  swap(maparray[0], maparray[1]);
1248 
1249  for(i = 3; i < P; i+=2)
1250  {
1251  signarray[i] = -1;
1252  }
1253  }
1254  }
1255  break;
1256  case 2:
1257  {
1258  for (i = 0; i < P; i++)
1259  {
1260  maparray[i] = order0+i;
1261  }
1262 
1263  if (edgeOrient == eBackwards)
1264  {
1265  swap(maparray[0], maparray[1]);
1266 
1267  for (i = 3; i < P; i+=2)
1268  {
1269  signarray[i] = -1;
1270  }
1271  }
1272  }
1273  break;
1274  case 3:
1275  {
1276  for (i = 0; i < P; i++)
1277  {
1278  maparray[i] = i*order0;
1279  }
1280 
1281  if (edgeOrient == eBackwards)
1282  {
1283  swap(maparray[0], maparray[1]);
1284 
1285  for (i = 3; i < P; i+=2)
1286  {
1287  signarray[i] = -1;
1288  }
1289  }
1290  }
1291  break;
1292  default:
1293  ASSERTL0(false, "eid must be between 0 and 3");
1294  break;
1295  }
1296  }
1297  else if(bType == LibUtilities::eGLL_Lagrange ||
1299  {
1300  switch (eid)
1301  {
1302  case 0:
1303  {
1304  for (i = 0; i < P; i++)
1305  {
1306  maparray[i] = i;
1307  }
1308  }
1309  break;
1310  case 1:
1311  {
1312  for (i = 0; i < P; i++)
1313  {
1314  maparray[i] = (i+1)*order0 - 1;
1315  }
1316  }
1317  break;
1318  case 2:
1319  {
1320  for (i = 0; i < P; i++)
1321  {
1322  maparray[i] = order0*(order1-1) + i;
1323  }
1324  }
1325  break;
1326  case 3:
1327  {
1328  for (i = 0; i < P; i++)
1329  {
1330  maparray[i] = order0*i;
1331  }
1332  }
1333  break;
1334  default:
1335  ASSERTL0(false, "eid must be between 0 and 3");
1336  break;
1337  }
1338  if (edgeOrient == eBackwards)
1339  {
1340  reverse(maparray.get(), maparray.get()+P);
1341  }
1342  }
1343  else
1344  {
1345  ASSERTL0(false, "Mapping not defined for this type of basis");
1346  }
1347 
1348  if (checkForZeroedModes)
1349  {
1350  if (bType == LibUtilities::eModified_A)
1351  {
1352  // Zero signmap and set maparray to zero if
1353  // elemental modes are not as large as face modesl
1354  for (int j = numModes; j < P; j++)
1355  {
1356  signarray[j] = 0.0;
1357  maparray[j] = maparray[0];
1358  }
1359  }
1360  else
1361  {
1362  ASSERTL0(false, "Different trace space edge dimension "
1363  "and element edge dimension not possible "
1364  "for GLL-Lagrange bases");
1365  }
1366  }
1367  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
Principle Modified Functions .
Definition: BasisType.h:48
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:55
Lagrange for SEM basis .
Definition: BasisType.h:54
Array< OneD, LibUtilities::BasisSharedPtr > m_base
LibUtilities::BasisType GetEdgeBasisType(const int i) const
This function returns the type of expansion basis on the i-th edge.
Definition: StdExpansion.h:412

◆ v_GetInteriorMap()

void Nektar::StdRegions::StdQuadExp::v_GetInteriorMap ( Array< OneD, unsigned int > &  outarray)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 853 of file StdQuadExp.cpp.

References ASSERTL0, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), Nektar::StdRegions::StdExpansion::m_base, and Nektar::StdRegions::StdExpansion::NumBndryCoeffs().

854  {
855  int i,j;
856  int cnt=0;
857  int nummodes0, nummodes1;
858  int startvalue = 0;
859  if(outarray.num_elements()!=GetNcoeffs()-NumBndryCoeffs())
860  {
861  outarray = Array<OneD, unsigned int>(GetNcoeffs()-NumBndryCoeffs());
862  }
863 
864  nummodes0 = m_base[0]->GetNumModes();
865  nummodes1 = m_base[1]->GetNumModes();
866 
867  const LibUtilities::BasisType Btype0 = GetBasisType(0);
868  const LibUtilities::BasisType Btype1 = GetBasisType(1);
869 
870  switch(Btype1)
871  {
873  startvalue = nummodes0;
874  break;
876  startvalue = 2*nummodes0;
877  break;
878  default:
879  ASSERTL0(0,"Mapping array is not defined for this expansion");
880  break;
881  }
882 
883  switch(Btype0)
884  {
886  startvalue++;
887  break;
889  startvalue+=2;
890  break;
891  default:
892  ASSERTL0(0,"Mapping array is not defined for this expansion");
893  break;
894  }
895 
896  for(i = 0; i < nummodes1-2; i++)
897  {
898  for(j = 0; j < nummodes0-2; j++)
899  {
900  outarray[cnt++]=startvalue+j;
901  }
902  startvalue+=nummodes0;
903  }
904  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
Principle Modified Functions .
Definition: BasisType.h:48
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:130
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:164
Lagrange for SEM basis .
Definition: BasisType.h:54
Array< OneD, LibUtilities::BasisSharedPtr > m_base

◆ v_GetNedges()

int Nektar::StdRegions::StdQuadExp::v_GetNedges ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 627 of file StdQuadExp.cpp.

628  {
629  return 4;
630  }

◆ v_GetNverts()

int Nektar::StdRegions::StdQuadExp::v_GetNverts ( ) const
protectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 622 of file StdQuadExp.cpp.

623  {
624  return 4;
625  }

◆ v_GetSimplexEquiSpacedConnectivity()

void Nektar::StdRegions::StdQuadExp::v_GetSimplexEquiSpacedConnectivity ( Array< OneD, int > &  conn,
bool  standard = true 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1811 of file StdQuadExp.cpp.

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

1814  {
1815  boost::ignore_unused(standard);
1816 
1817  int np1 = m_base[0]->GetNumPoints();
1818  int np2 = m_base[1]->GetNumPoints();
1819  int np = max(np1,np2);
1820 
1821  conn = Array<OneD, int>(6*(np-1)*(np-1));
1822 
1823  int row = 0;
1824  int rowp1 = 0;
1825  int cnt = 0;
1826  for(int i = 0; i < np-1; ++i)
1827  {
1828  rowp1 += np;
1829  for(int j = 0; j < np-1; ++j)
1830  {
1831  conn[cnt++] = row +j;
1832  conn[cnt++] = row +j+1;
1833  conn[cnt++] = rowp1 +j;
1834 
1835  conn[cnt++] = rowp1 +j+1;
1836  conn[cnt++] = rowp1 +j;
1837  conn[cnt++] = row +j+1;
1838  }
1839  row += np;
1840  }
1841  }
Array< OneD, LibUtilities::BasisSharedPtr > m_base

◆ v_GetVertexMap()

int Nektar::StdRegions::StdQuadExp::v_GetVertexMap ( int  localVertexId,
bool  useCoeffPacking = false 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 906 of file StdQuadExp.cpp.

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

907  {
908  int localDOF = 0;
909 
910  if(useCoeffPacking == true)
911  {
912  switch(localVertexId)
913  {
914  case 0:
915  {
916  localDOF = 0;
917  }
918  break;
919  case 1:
920  {
922  {
923  localDOF = m_base[0]->GetNumModes()-1;
924  }
925  else
926  {
927  localDOF = 1;
928  }
929  }
930  break;
931  case 2:
932  {
934  {
935  localDOF = m_base[0]->GetNumModes() * (m_base[1]->GetNumModes()-1);
936  }
937  else
938  {
939  localDOF = m_base[0]->GetNumModes();
940  }
941  }
942  break;
943  case 3:
944  {
946  {
947  localDOF = m_base[0]->GetNumModes()*m_base[1]->GetNumModes()-1;
948  }
949  else
950  {
951  localDOF = m_base[0]->GetNumModes()+1;
952  }
953  }
954  break;
955  default:
956  ASSERTL0(false,"eid must be between 0 and 3");
957  break;
958  }
959 
960  }
961  else
962  {
963  switch(localVertexId)
964  {
965  case 0:
966  {
967  localDOF = 0;
968  }
969  break;
970  case 1:
971  {
973  {
974  localDOF = m_base[0]->GetNumModes()-1;
975  }
976  else
977  {
978  localDOF = 1;
979  }
980  }
981  break;
982  case 2:
983  {
985  {
986  localDOF = m_base[0]->GetNumModes()*m_base[1]->GetNumModes()-1;
987  }
988  else
989  {
990  localDOF = m_base[0]->GetNumModes()+1;
991  }
992  }
993  break;
994  case 3:
995  {
997  {
998  localDOF = m_base[0]->GetNumModes() * (m_base[1]->GetNumModes()-1);
999  }
1000  else
1001  {
1002  localDOF = m_base[0]->GetNumModes();
1003  }
1004  }
1005  break;
1006  default:
1007  ASSERTL0(false,"eid must be between 0 and 3");
1008  break;
1009  }
1010  }
1011  return localDOF;
1012  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:164
Lagrange for SEM basis .
Definition: BasisType.h:54
Array< OneD, LibUtilities::BasisSharedPtr > m_base

◆ v_HelmholtzMatrixOp()

void Nektar::StdRegions::StdQuadExp::v_HelmholtzMatrixOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 1778 of file StdQuadExp.cpp.

References Nektar::StdRegions::StdExpansion2D::v_HelmholtzMatrixOp_MatFree().

1781  {
1782  StdQuadExp::v_HelmholtzMatrixOp_MatFree(inarray,outarray,mkey);
1783  }
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)

◆ v_Integral()

NekDouble Nektar::StdRegions::StdQuadExp::v_Integral ( const Array< OneD, const NekDouble > &  inarray)
protectedvirtual

Integrates the specified function over the domain.

See also
StdRegions::StdExpansion::Integral.

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 80 of file StdQuadExp.cpp.

References Nektar::StdRegions::StdExpansion2D::Integral(), and Nektar::StdRegions::StdExpansion::m_base.

82  {
83  Array<OneD, const NekDouble> w0 = m_base[0]->GetW();
84  Array<OneD, const NekDouble> w1 = m_base[1]->GetW();
85 
86  return StdExpansion2D::Integral(inarray,w0,w1);
87  }
NekDouble Integral(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &w0, const Array< OneD, const NekDouble > &w1)
Array< OneD, LibUtilities::BasisSharedPtr > m_base

◆ v_IProductWRTBase()

void Nektar::StdRegions::StdQuadExp::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 and put into outarray.

\( \begin{array}{rcl} I_{pq} = (\phi_q \phi_q, u) & = & \sum_{i=0}^{nq_0} \sum_{j=0}^{nq_1} \phi_p(\xi_{0,i}) \phi_q(\xi_{1,j}) w^0_i w^1_j u(\xi_{0,i} \xi_{1,j}) \\ & = & \sum_{i=0}^{nq_0} \phi_p(\xi_{0,i}) \sum_{j=0}^{nq_1} \phi_q(\xi_{1,j}) \tilde{u}_{i,j} \end{array} \)

where

\( \tilde{u}_{i,j} = w^0_i w^1_j u(\xi_{0,i},\xi_{1,j}) \)

which can be implemented as

\( f_{qi} = \sum_{j=0}^{nq_1} \phi_q(\xi_{1,j}) \tilde{u}_{i,j} = {\bf B_1 U} \) \( I_{pq} = \sum_{i=0}^{nq_0} \phi_p(\xi_{0,i}) f_{qi} = {\bf B_0 F} \)

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 385 of file StdQuadExp.cpp.

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

Referenced by v_FwdTrans().

388  {
389  if(m_base[0]->Collocation() && m_base[1]->Collocation())
390  {
391  MultiplyByQuadratureMetric(inarray,outarray);
392  }
393  else
394  {
395  StdQuadExp::v_IProductWRTBase_SumFac(inarray,outarray);
396  }
397  }
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:945
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
Definition: StdQuadExp.cpp:399
Array< OneD, LibUtilities::BasisSharedPtr > m_base

◆ v_IProductWRTBase_MatOp()

void Nektar::StdRegions::StdQuadExp::v_IProductWRTBase_MatOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 428 of file StdQuadExp.cpp.

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

431  {
432  int nq = GetTotPoints();
433  StdMatrixKey iprodmatkey(eIProductWRTBase,DetShapeType(),*this);
434  DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
435 
436  Blas::Dgemv('N',m_ncoeffs,nq,1.0,iprodmat->GetPtr().get(),
437  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
438  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:469
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:714
static void Dgemv(const char &trans, const int &m, const int &n, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = A x where A[m x n].
Definition: Blas.hpp:168
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140

◆ v_IProductWRTBase_SumFac()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 399 of file StdQuadExp.cpp.

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

Referenced by v_IProductWRTBase().

403  {
404  int nquad0 = m_base[0]->GetNumPoints();
405  int nquad1 = m_base[1]->GetNumPoints();
406  int order0 = m_base[0]->GetNumModes();
407 
408  if(multiplybyweights)
409  {
410  Array<OneD,NekDouble> tmp(nquad0*nquad1+nquad1*order0);
411  Array<OneD,NekDouble> wsp(tmp+nquad0*nquad1);
412 
413  // multiply by integration constants
414  MultiplyByQuadratureMetric(inarray,tmp);
415  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
416  m_base[1]->GetBdata(),
417  tmp,outarray,wsp,true,true);
418  }
419  else
420  {
421  Array<OneD,NekDouble> wsp(nquad1*order0);
422  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
423  m_base[1]->GetBdata(),
424  inarray,outarray,wsp,true,true);
425  }
426  }
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:945
void IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
Array< OneD, LibUtilities::BasisSharedPtr > m_base

◆ v_IProductWRTBase_SumFacKernel()

void Nektar::StdRegions::StdQuadExp::v_IProductWRTBase_SumFacKernel ( const Array< OneD, const NekDouble > &  base0,
const Array< OneD, const NekDouble > &  base1,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
Array< OneD, NekDouble > &  wsp,
bool  doCheckCollDir0,
bool  doCheckCollDir1 
)
protectedvirtual

Implements Nektar::StdRegions::StdExpansion2D.

Definition at line 513 of file StdQuadExp.cpp.

References ASSERTL1, Blas::Ddot(), Blas::Dgemm(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, and Vmath::Vcopy().

521  {
522  int nquad0 = m_base[0]->GetNumPoints();
523  int nquad1 = m_base[1]->GetNumPoints();
524  int nmodes0 = m_base[0]->GetNumModes();
525  int nmodes1 = m_base[1]->GetNumModes();
526 
527  bool colldir0 = doCheckCollDir0?(m_base[0]->Collocation()):false;
528  bool colldir1 = doCheckCollDir1?(m_base[1]->Collocation()):false;
529 
530  if(colldir0 && colldir1)
531  {
532  Vmath::Vcopy(m_ncoeffs,inarray.get(),1,outarray.get(),1);
533  }
534  else if(colldir0)
535  {
536  Blas::Dgemm('N', 'N',nmodes0,nmodes1, nquad1,1.0, inarray.get(),
537  nmodes0, base1.get(), nquad1, 0.0,outarray.get(),nmodes0);
538  }
539  else if(colldir1)
540  {
541  Blas::Dgemm('T','N',nmodes0,nquad1,nquad0,1.0,base0.get(),
542  nquad0,inarray.get(),nquad0,0.0,outarray.get(),nmodes0);
543  }
544  else
545  {
546  ASSERTL1(wsp.num_elements()>=nquad1*nmodes0,"Workspace size is not sufficient");
547 
548 #if 1
549  Blas::Dgemm('T','N',nmodes0,nquad1,nquad0,1.0,base0.get(),
550  nquad0,inarray.get(),nquad0,0.0,wsp.get(),nmodes0);
551 
552 #else
553  for(int i = 0; i < nmodes0; ++i)
554  {
555  for(int j = 0; j < nquad1; ++j)
556  {
557  wsp[j*nmodes0+i] = Blas::Ddot(nquad0,
558  base0.get()+i*nquad0,1,
559  inarray.get()+j*nquad0,1);
560  }
561  }
562 #endif
563  Blas::Dgemm('N', 'N',nmodes0,nmodes1, nquad1,1.0, wsp.get(),
564  nmodes0, base1.get(), nquad1, 0.0,outarray.get(),nmodes0);
565  }
566  }
static void Dgemm(const char &transa, const char &transb, const int &m, const int &n, const int &k, const double &alpha, const double *a, const int &lda, const double *b, const int &ldb, const double &beta, double *c, const int &ldc)
BLAS level 3: Matrix-matrix multiply C = A x B where A[m x n], B[n x k], C[m x k].
Definition: Blas.hpp:213
static double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = .
Definition: Blas.hpp:140
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064

◆ v_IProductWRTDerivBase()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 440 of file StdQuadExp.cpp.

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

443  {
444  StdQuadExp::IProductWRTDerivBase_SumFac(dir,inarray,outarray);
445  }
void IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)

◆ v_IProductWRTDerivBase_MatOp()

void Nektar::StdRegions::StdQuadExp::v_IProductWRTDerivBase_MatOp ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 478 of file StdQuadExp.cpp.

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

481  {
482  ASSERTL0((dir==0)||(dir==1),"input dir is out of range");
483 
484  int nq = GetTotPoints();
485  MatrixType mtype;
486 
487  if(dir) // dir == 1
488  {
489  mtype = eIProductWRTDerivBase1;
490  }
491  else // dir == 0
492  {
493  mtype = eIProductWRTDerivBase0;
494  }
495 
496  StdMatrixKey iprodmatkey(mtype,DetShapeType(),*this);
497  DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
498 
499  Blas::Dgemv('N',m_ncoeffs,nq,1.0,iprodmat->GetPtr().get(),
500  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
501  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:469
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:714
static void Dgemv(const char &trans, const int &m, const int &n, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = A x where A[m x n].
Definition: Blas.hpp:168
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140

◆ v_IProductWRTDerivBase_SumFac()

void Nektar::StdRegions::StdQuadExp::v_IProductWRTDerivBase_SumFac ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 447 of file StdQuadExp.cpp.

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

450  {
451  ASSERTL0((dir==0)||(dir==1),"input dir is out of range");
452 
453  int nquad0 = m_base[0]->GetNumPoints();
454  int nquad1 = m_base[1]->GetNumPoints();
455  int nqtot = nquad0*nquad1;
456  int order0 = m_base[0]->GetNumModes();
457 
458  Array<OneD,NekDouble> tmp(nqtot+nquad1*order0);
459  Array<OneD,NekDouble> wsp(tmp+nqtot);
460 
461  // multiply by integration constants
462  MultiplyByQuadratureMetric(inarray,tmp);
463 
464  if(dir) // dir == 1
465  {
466  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
467  m_base[1]->GetDbdata(),
468  tmp,outarray,wsp,true,false);
469  }
470  else // dir == 0
471  {
472  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
473  m_base[1]->GetBdata(),
474  tmp,outarray,wsp,false,true);
475  }
476  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:945
void IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
Array< OneD, LibUtilities::BasisSharedPtr > m_base

◆ v_IsBoundaryInteriorExpansion()

bool Nektar::StdRegions::StdQuadExp::v_IsBoundaryInteriorExpansion ( )
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 747 of file StdQuadExp.cpp.

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

748  {
749  bool returnval = false;
750 
753  {
756  {
757  returnval = true;
758  }
759  }
760 
761  return returnval;
762  }
Principle Modified Functions .
Definition: BasisType.h:48
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:164
Lagrange for SEM basis .
Definition: BasisType.h:54
Array< OneD, LibUtilities::BasisSharedPtr > m_base

◆ v_LaplacianMatrixOp() [1/2]

void Nektar::StdRegions::StdQuadExp::v_LaplacianMatrixOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 1754 of file StdQuadExp.cpp.

References Nektar::StdRegions::StdExpansion2D::v_LaplacianMatrixOp_MatFree().

1758  {
1759  StdQuadExp::v_LaplacianMatrixOp_MatFree(inarray,outarray,mkey);
1760  }
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)

◆ v_LaplacianMatrixOp() [2/2]

void Nektar::StdRegions::StdQuadExp::v_LaplacianMatrixOp ( const int  k1,
const int  k2,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 1762 of file StdQuadExp.cpp.

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

1766  {
1767  StdExpansion::LaplacianMatrixOp_MatFree(k1,k2,inarray,outarray,mkey);
1768  }
void LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)

◆ v_LocCoordToLocCollapsed()

void Nektar::StdRegions::StdQuadExp::v_LocCoordToLocCollapsed ( const Array< OneD, const NekDouble > &  xi,
Array< OneD, NekDouble > &  eta 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 573 of file StdQuadExp.cpp.

575  {
576  eta[0] = xi[0];
577  eta[1] = xi[1];
578  }

◆ v_MassMatrixOp()

void Nektar::StdRegions::StdQuadExp::v_MassMatrixOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 1746 of file StdQuadExp.cpp.

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

1750  {
1751  StdExpansion::MassMatrixOp_MatFree(inarray,outarray,mkey);
1752  }
void MassMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)

◆ v_MultiplyByStdQuadratureMetric()

void Nektar::StdRegions::StdQuadExp::v_MultiplyByStdQuadratureMetric ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1786 of file StdQuadExp.cpp.

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

1789  {
1790  int i;
1791  int nquad0 = m_base[0]->GetNumPoints();
1792  int nquad1 = m_base[1]->GetNumPoints();
1793 
1794  const Array<OneD, const NekDouble>& w0 = m_base[0]->GetW();
1795  const Array<OneD, const NekDouble>& w1 = m_base[1]->GetW();
1796 
1797  // multiply by integration constants
1798  for(i = 0; i < nquad1; ++i)
1799  {
1800  Vmath::Vmul(nquad0, inarray.get()+i*nquad0,1,
1801  w0.get(),1,outarray.get()+i*nquad0,1);
1802  }
1803 
1804  for(i = 0; i < nquad0; ++i)
1805  {
1806  Vmath::Vmul(nquad1,outarray.get()+i,nquad0,w1.get(),1,
1807  outarray.get()+i,nquad0);
1808  }
1809  }
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:186

◆ v_NumBndryCoeffs()

int Nektar::StdRegions::StdQuadExp::v_NumBndryCoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 709 of file StdQuadExp.cpp.

References ASSERTL1, Nektar::LibUtilities::eGauss_Lagrange, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisNumModes(), and Nektar::StdRegions::StdExpansion::GetBasisType().

710  {
714  "BasisType is not a boundary interior form");
718  "BasisType is not a boundary interior form");
719 
720  return 4 + 2*(GetBasisNumModes(0)-2) + 2*(GetBasisNumModes(1)-2);
721  }
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:177
Principle Modified Functions .
Definition: BasisType.h:48
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:55
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:164
Lagrange for SEM basis .
Definition: BasisType.h:54
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250

◆ v_NumDGBndryCoeffs()

int Nektar::StdRegions::StdQuadExp::v_NumDGBndryCoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 723 of file StdQuadExp.cpp.

References ASSERTL1, Nektar::LibUtilities::eGauss_Lagrange, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisNumModes(), and Nektar::StdRegions::StdExpansion::GetBasisType().

724  {
728  "BasisType is not a boundary interior form");
732  "BasisType is not a boundary interior form");
733 
734  return 2*GetBasisNumModes(0) + 2*GetBasisNumModes(1);
735  }
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:177
Principle Modified Functions .
Definition: BasisType.h:48
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:55
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:164
Lagrange for SEM basis .
Definition: BasisType.h:54
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250

◆ v_PhysDeriv() [1/2]

void Nektar::StdRegions::StdQuadExp::v_PhysDeriv ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d0,
Array< OneD, NekDouble > &  out_d1,
Array< OneD, NekDouble > &  out_d2 = NullNekDouble1DArray 
)
protectedvirtual

Calculate the derivative of the physical points.

For quadrilateral region can use the Tensor_Deriv function defined under StdExpansion.

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 99 of file StdQuadExp.cpp.

References Nektar::StdRegions::StdExpansion2D::PhysTensorDeriv().

Referenced by v_StdPhysDeriv().

103  {
104  boost::ignore_unused(out_d2);
105  PhysTensorDeriv(inarray, out_d0, out_d1);
106  }
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d0, Array< OneD, NekDouble > &outarray_d1)
Calculate the 2D derivative in the local tensor/collapsed coordinate at the physical points...

◆ v_PhysDeriv() [2/2]

void Nektar::StdRegions::StdQuadExp::v_PhysDeriv ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d0 
)
protectedvirtual

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

See also
StdRegions::StdExpansion::PhysDeriv

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 108 of file StdQuadExp.cpp.

References ASSERTL1, Nektar::NullNekDouble1DArray, and Nektar::StdRegions::StdExpansion2D::PhysTensorDeriv().

111  {
112  switch(dir)
113  {
114  case 0:
115  {
116  PhysTensorDeriv(inarray, outarray, NullNekDouble1DArray);
117  }
118  break;
119  case 1:
120  {
121  PhysTensorDeriv(inarray, NullNekDouble1DArray, outarray);
122  }
123  break;
124  default:
125  {
126  ASSERTL1(false,"input dir is out of range");
127  }
128  break;
129  }
130  }
static Array< OneD, NekDouble > NullNekDouble1DArray
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d0, Array< OneD, NekDouble > &outarray_d1)
Calculate the 2D derivative in the local tensor/collapsed coordinate at the physical points...
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250

◆ v_ReduceOrderCoeffs()

void Nektar::StdRegions::StdQuadExp::v_ReduceOrderCoeffs ( int  numMin,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 1696 of file StdQuadExp.cpp.

References Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eOrtho_A, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::LibUtilities::InterpCoeff2D(), Nektar::StdRegions::StdExpansion::m_base, Vmath::Vcopy(), and Vmath::Zero().

1700  {
1701  int n_coeffs = inarray.num_elements();
1702 
1703 
1704  Array<OneD, NekDouble> coeff(n_coeffs);
1705  Array<OneD, NekDouble> coeff_tmp(n_coeffs,0.0);
1706  Array<OneD, NekDouble> tmp;
1707  Array<OneD, NekDouble> tmp2;
1708 
1709  int nmodes0 = m_base[0]->GetNumModes();
1710  int nmodes1 = m_base[1]->GetNumModes();
1711  int numMax = nmodes0;
1712 
1713  Vmath::Vcopy(n_coeffs,inarray,1,coeff_tmp,1);
1714 
1715  const LibUtilities::PointsKey Pkey0(
1717  const LibUtilities::PointsKey Pkey1(
1719 
1720  LibUtilities::BasisKey b0(m_base[0]->GetBasisType(),nmodes0,Pkey0);
1721  LibUtilities::BasisKey b1(m_base[1]->GetBasisType(),nmodes1,Pkey1);
1722 
1723  LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A,nmodes0,Pkey0);
1724  LibUtilities::BasisKey bortho1(LibUtilities::eOrtho_A,nmodes1,Pkey1);
1725 
1727  b0, b1, coeff_tmp, bortho0, bortho1, coeff);
1728 
1729  Vmath::Zero(n_coeffs,coeff_tmp,1);
1730 
1731  int cnt = 0;
1732  for (int i = 0; i < numMin+1; ++i)
1733  {
1734  Vmath::Vcopy(numMin,
1735  tmp = coeff+cnt,1,
1736  tmp2 = coeff_tmp+cnt,1);
1737 
1738  cnt = i*numMax;
1739  }
1740 
1742  bortho0, bortho1, coeff_tmp, b0, b1, outarray);
1743  }
Principle Orthogonal Functions .
Definition: BasisType.h:45
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:164
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:376
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void InterpCoeff2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
Definition: InterpCoeff.cpp:72
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51

◆ v_StdPhysDeriv() [1/2]

void Nektar::StdRegions::StdQuadExp::v_StdPhysDeriv ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d0,
Array< OneD, NekDouble > &  out_d1,
Array< OneD, NekDouble > &  out_d2 = NullNekDouble1DArray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 132 of file StdQuadExp.cpp.

References v_PhysDeriv().

136  {
137  boost::ignore_unused(out_d2);
138  //PhysTensorDeriv(inarray, out_d0, out_d1);
139  StdQuadExp::v_PhysDeriv(inarray, out_d0, out_d1);
140  }
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=NullNekDouble1DArray)
Calculate the derivative of the physical points.
Definition: StdQuadExp.cpp:99

◆ v_StdPhysDeriv() [2/2]

void Nektar::StdRegions::StdQuadExp::v_StdPhysDeriv ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 142 of file StdQuadExp.cpp.

References v_PhysDeriv().

145  {
146  //PhysTensorDeriv(inarray, outarray);
147  StdQuadExp::v_PhysDeriv(dir,inarray,outarray);
148 
149  }
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=NullNekDouble1DArray)
Calculate the derivative of the physical points.
Definition: StdQuadExp.cpp:99

◆ v_SVVLaplacianFilter()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 1544 of file StdQuadExp.cpp.

References Nektar::StdRegions::StdExpansion::BwdTrans(), Nektar::StdRegions::StdMatrixKey::ConstFactorExists(), Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::eFactorSVVDGKerDiffCoeff, Nektar::StdRegions::eFactorSVVDiffCoeff, Nektar::StdRegions::eFactorSVVPowerKerDiffCoeff, Nektar::LibUtilities::eOrtho_A, Nektar::StdRegions::StdExpansion::FwdTrans(), Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::StdRegions::kSVVDGFilter, Nektar::StdRegions::kSVVDGFiltermodesmax, Nektar::StdRegions::kSVVDGFiltermodesmin, and Nektar::StdRegions::StdExpansion::m_base.

1546  {
1547  // Generate an orthonogal expansion
1548  int qa = m_base[0]->GetNumPoints();
1549  int qb = m_base[1]->GetNumPoints();
1550  int nmodes_a = m_base[0]->GetNumModes();
1551  int nmodes_b = m_base[1]->GetNumModes();
1552  int nmodes = min(nmodes_a,nmodes_b);
1553  // Declare orthogonal basis.
1554  LibUtilities::PointsKey pa(qa,m_base[0]->GetPointsType());
1555  LibUtilities::PointsKey pb(qb,m_base[1]->GetPointsType());
1556 
1557  LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A,nmodes_a,pa);
1558  LibUtilities::BasisKey Bb(LibUtilities::eOrtho_A,nmodes_b,pb);
1559  StdQuadExp OrthoExp(Ba,Bb);
1560 
1561  //SVV parameters loaded from the .xml case file
1562  Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
1563 
1564  // project onto modal space.
1565  OrthoExp.FwdTrans(array,orthocoeffs);
1566 
1567  if(mkey.ConstFactorExists(eFactorSVVPowerKerDiffCoeff)) // Rodrigo's power kernel
1568  {
1569  NekDouble cutoff = mkey.GetConstFactor(eFactorSVVCutoffRatio);
1570  NekDouble SvvDiffCoeff =
1571  mkey.GetConstFactor(eFactorSVVPowerKerDiffCoeff)*
1572  mkey.GetConstFactor(eFactorSVVDiffCoeff);
1573 
1574  for(int j = 0; j < nmodes_a; ++j)
1575  {
1576  for(int k = 0; k < nmodes_b; ++k)
1577  {
1578  // linear space but makes high modes very negative
1579  NekDouble fac = std::max(
1580  pow((1.0*j)/(nmodes_a-1),cutoff*nmodes_a),
1581  pow((1.0*k)/(nmodes_b-1),cutoff*nmodes_b));
1582 
1583  orthocoeffs[j*nmodes_b+k] *= SvvDiffCoeff * fac;
1584  }
1585  }
1586  }
1587  else if(mkey.ConstFactorExists(eFactorSVVDGKerDiffCoeff)) // Rodrigo/mansoor's DG kernel
1588  {
1589  NekDouble SvvDiffCoeff =
1590  mkey.GetConstFactor(eFactorSVVDGKerDiffCoeff)*
1591  mkey.GetConstFactor(eFactorSVVDiffCoeff);
1592  int max_ab = max(nmodes_a-kSVVDGFiltermodesmin,
1593  nmodes_b-kSVVDGFiltermodesmin);
1594  max_ab = max(max_ab,0);
1595  max_ab = min(max_ab,kSVVDGFiltermodesmax-kSVVDGFiltermodesmin);
1596 
1597  for(int j = 0; j < nmodes_a; ++j)
1598  {
1599  for(int k = 0; k < nmodes_b; ++k)
1600  {
1601  int maxjk = max(j,k);
1602  maxjk = min(maxjk,kSVVDGFiltermodesmax-1);
1603 
1604  orthocoeffs[j*nmodes_b+k] *= SvvDiffCoeff *
1605  kSVVDGFilter[max_ab][maxjk];
1606  }
1607  }
1608  }
1609  else
1610  {
1611  NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDiffCoeff);
1612  //Exponential Kernel implementation
1613  int cutoff = (int) (mkey.GetConstFactor(eFactorSVVCutoffRatio)*
1614  min(nmodes_a,nmodes_b));
1615 
1616  //counters for scanning through orthocoeffs array
1617  int cnt = 0;
1618 
1619  //------"New" Version August 22nd '13--------------------
1620  for(int j = 0; j < nmodes_a; ++j)
1621  {
1622  for(int k = 0; k < nmodes_b; ++k)
1623  {
1624  if(j + k >= cutoff)//to filter out only the "high-modes"
1625  {
1626  orthocoeffs[j*nmodes_b+k] *=
1627  (SvvDiffCoeff*exp(-(j+k-nmodes)*(j+k-nmodes)/
1628  ((NekDouble)((j+k-cutoff+1)*
1629  (j+k-cutoff+1)))));
1630  }
1631  else
1632  {
1633  orthocoeffs[j*nmodes_b+k] *= 0.0;
1634  }
1635  cnt++;
1636  }
1637  }
1638  }
1639 
1640  // backward transform to physical space
1641  OrthoExp.BwdTrans(orthocoeffs,array);
1642  }
const int kSVVDGFiltermodesmax
Definition: StdRegions.hpp:385
const int kSVVDGFiltermodesmin
Definition: StdRegions.hpp:384
Principle Orthogonal Functions .
Definition: BasisType.h:45
double NekDouble
const NekDouble kSVVDGFilter[9][11]
Definition: StdRegions.hpp:387
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:215
Array< OneD, LibUtilities::BasisSharedPtr > m_base

◆ v_WeakDerivMatrixOp()

void Nektar::StdRegions::StdQuadExp::v_WeakDerivMatrixOp ( const int  i,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 1770 of file StdQuadExp.cpp.

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

1774  {
1775  StdExpansion::WeakDerivMatrixOp_MatFree(i,inarray,outarray,mkey);
1776  }
void WeakDerivMatrixOp_MatFree(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)