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 GetTraceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th trace. More...
 
int GetTraceIntNcoeffs (const int i) const
 
int GetTraceNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th trace. More...
 
const LibUtilities::BasisKey GetTraceBasisKey (const int i, int k=-1) const
 This function returns the basis key belonging to the i-th trace. More...
 
LibUtilities::PointsKey GetTracePointsKey (const int i, int k=-1) const
 This function returns the basis key belonging to the i-th trace. More...
 
int NumBndryCoeffs (void) const
 
int NumDGBndryCoeffs (void) const
 
const LibUtilities::PointsKey GetNodalPointsKey () const
 This function returns the type of expansion Nodal point type if defined. More...
 
int GetNtraces () const
 Returns the number of trace elements connected to this element. More...
 
LibUtilities::ShapeType DetShapeType () const
 This function returns the shape of the expansion domain. More...
 
std::shared_ptr< StdExpansionGetStdExp (void) const
 
std::shared_ptr< StdExpansionGetLinStdExp (void) const
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion ()
 
bool IsNodalNonTensorialExp ()
 
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space. More...
 
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space. More...
 
void FwdTrans_BndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble Integral (const Array< OneD, const NekDouble > &inarray)
 This function integrates the specified function over the domain. More...
 
void FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 This function fills the array outarray with the mode-th mode of the expansion. More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 this function calculates the inner product of a given function f with the different modes of the expansion More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check)
 
void IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
int GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void SetElmtId (const int id)
 Set the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2=NullNekDouble1DArray, Array< OneD, NekDouble > &coords_3=NullNekDouble1DArray)
 this function returns the physical coordinates of the quadrature points of the expansion More...
 
void GetCoord (const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
 given the coordinates of a point of the element in the local collapsed coordinate system, this function calculates the physical coordinates of the point More...
 
DNekMatSharedPtr GetStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr GetStdStaticCondMatrix (const StdMatrixKey &mkey)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
void DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
int CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
NekDouble StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
int GetCoordim ()
 
void GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
void GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
int GetVertexMap (const int localVertexId, bool useCoeffPacking=false)
 
void GetTraceToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
 
void GetTraceInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eForwards)
 
void GetTraceNumModes (const int tid, int &numModes0, int &numModes1, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2)
 
void MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
DNekMatSharedPtr CreateGeneralMatrix (const StdMatrixKey &mkey)
 this function generates the mass matrix \(\mathbf{M}[i][j] = \int \phi_i(\mathbf{x}) \phi_j(\mathbf{x}) d\mathbf{x}\) More...
 
void GeneralMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
void ExponentialFilter (Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
 
void LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
DNekMatSharedPtr GenMatrix (const StdMatrixKey &mkey)
 
void PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void PhysDeriv_s (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds)
 
void PhysDeriv_n (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn)
 
void PhysDirectionalDeriv (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &outarray)
 
void StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluateBasis (const Array< OneD, const NekDouble > &coords, int mode)
 This function evaluates the basis function mode mode at a point coords of the domain. More...
 
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta. More...
 
void LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi)
 Convert local collapsed coordinates eta into local cartesian coordinate xi. More...
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
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)
 
NekDouble Linf (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( L_\infty\) error \( |\epsilon|_\infty = \max |u - u_{exact}|\) where \( u_{exact}\) is given by the array sol. More...
 
NekDouble L2 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( L_2\) error, \( | \epsilon |_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 dx \right]^{1/2} d\xi_1 \) where \( u_{exact}\) is given by the array sol. More...
 
NekDouble H1 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( H^1\) error, \( | \epsilon |^1_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 + \nabla(u - u_{exact})\cdot\nabla(u - u_{exact})\cdot dx \right]^{1/2} d\xi_1 \) where \( u_{exact}\) is given by the array sol. More...
 
const LibUtilities::PointsKeyVector GetPointsKeys () const
 
DNekMatSharedPtr BuildInverseTransformationMatrix (const DNekScalMatSharedPtr &m_transformationmatrix)
 
void PhysInterpToSimplexEquiSpaced (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int npset=-1)
 This function performs an interpolation from the physical space points provided at input into an array of equispaced points which are not the collapsed coordinate. So for a tetrahedron you will only get a tetrahedral number of values. More...
 
void GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true)
 This function provides the connectivity of local simplices (triangles or tets) to connect the equispaced data points provided by PhysInterpToSimplexEquiSpaced. More...
 
void EquiSpacedToCoeffs (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs a projection/interpolation from the equispaced points sometimes used in post-processing onto the coefficient space. More...
 
template<class T >
std::shared_ptr< T > as ()
 
void IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
void GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 

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_LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi)
 
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_GetNtraces () const
 
virtual int v_GetTraceNcoeffs (const int i) const
 
virtual int v_GetTraceNumPoints (const int i) const
 
virtual int v_NumBndryCoeffs () const
 
virtual int v_NumDGBndryCoeffs () const
 
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
virtual const LibUtilities::BasisKey v_GetTraceBasisKey (const int i, const int j) 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 NekDouble v_PhysEvaluateBasis (const Array< OneD, const NekDouble > &coords, int mode)
 This function evaluates the basis function mode mode at a point coords of the domain. More...
 
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_GetTraceToElementMap (const int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient=eForwards, int P=-1, int Q=-1)
 
virtual void v_GetTraceInteriorToElementMap (const int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation edgeOrient=eForwards)
 
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 void v_GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion
DNekMatSharedPtr CreateStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr CreateStdStaticCondMatrix (const StdMatrixKey &mkey)
 Create the static condensation of a matrix when using a boundary interior decomposition. More...
 
void BwdTrans_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase_SumFac (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GeneralMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
 
void LaplacianMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp_MatFree (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void HelmholtzMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_SetCoeffsToOrientation (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)
 
template<int DIR>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals)
 This function performs the barycentric interpolation of the polynomial stored in coord at a point physvals using barycentric interpolation weights in direction. More...
 
template<int DIR>
NekDouble BaryEvaluateBasis (const NekDouble &coord, const int &mode)
 

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
 

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  }

◆ ~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.

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:221
Array< OneD, LibUtilities::BasisSharedPtr > m_base
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdQuadExp.cpp:170
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199

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

◆ 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.

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)

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

Referenced by v_BwdTrans().

◆ 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.

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.size()>=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  }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
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 op(A)[m x k], op(B)[k x n], C[m x n] DGEMM perfo...
Definition: Blas.hpp:394

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

◆ 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 717 of file StdQuadExp.cpp.

720  {
721  int nmodes = nummodes[modes_offset]*nummodes[modes_offset+1];
722  modes_offset += 2;
723 
724  return nmodes;
725  }

◆ 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 1472 of file StdQuadExp.cpp.

1473  {
1474  return GenMatrix(mkey);
1475  }
DNekMatSharedPtr GenMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:850

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

◆ v_DetShapeType()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 683 of file StdQuadExp.cpp.

684  {
686  }

References Nektar::LibUtilities::eQuadrilateral.

◆ 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 1605 of file StdQuadExp.cpp.

1610  {
1611  // Generate an orthogonal expansion
1612  int qa = m_base[0]->GetNumPoints();
1613  int qb = m_base[1]->GetNumPoints();
1614  int nmodesA = m_base[0]->GetNumModes();
1615  int nmodesB = m_base[1]->GetNumModes();
1616  int P = nmodesA - 1;
1617  int Q = nmodesB - 1;
1618 
1619  // Declare orthogonal basis.
1620  LibUtilities::PointsKey pa(qa,m_base[0]->GetPointsType());
1621  LibUtilities::PointsKey pb(qb,m_base[1]->GetPointsType());
1622 
1623  LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A, nmodesA, pa);
1624  LibUtilities::BasisKey Bb(LibUtilities::eOrtho_A, nmodesB, pb);
1625  StdQuadExp OrthoExp(Ba,Bb);
1626 
1627  // Cutoff
1628  int Pcut = cutoff*P;
1629  int Qcut = cutoff*Q;
1630 
1631  // Project onto orthogonal space.
1632  Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
1633  OrthoExp.FwdTrans(array,orthocoeffs);
1634 
1635  //
1636  NekDouble fac, fac1, fac2;
1637  for(int i = 0; i < nmodesA; ++i)
1638  {
1639  for(int j = 0; j < nmodesB; ++j)
1640  {
1641  //to filter out only the "high-modes"
1642  if(i > Pcut || j > Qcut)
1643  {
1644  fac1 = (NekDouble) (i - Pcut)/( (NekDouble)(P - Pcut) );
1645  fac2 = (NekDouble) (j - Qcut)/( (NekDouble)(Q - Qcut) );
1646  fac = max(fac1, fac2);
1647  fac = pow(fac, exponent);
1648  orthocoeffs[i*nmodesB+j] *= exp(-alpha*fac);
1649  }
1650  }
1651  }
1652 
1653  // backward transform to physical space
1654  OrthoExp.BwdTrans(orthocoeffs,array);
1655  }
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:208
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:45
double NekDouble
P
Definition: main.py:133

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 main::P.

◆ 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 593 of file StdQuadExp.cpp.

595  {
596  int i;
597  int nquad0 = m_base[0]->GetNumPoints();
598  int nquad1 = m_base[1]->GetNumPoints();
599  Array<OneD, const NekDouble> base0 = m_base[0]->GetBdata();
600  Array<OneD, const NekDouble> base1 = m_base[1]->GetBdata();
601  int btmp0 = m_base[0]->GetNumModes();
602  int mode0 = mode%btmp0;
603  int mode1 = mode/btmp0;
604 
605  ASSERTL2(mode1 == (int)floor((1.0*mode)/btmp0),
606  "Integer Truncation not Equiv to Floor");
607 
608  ASSERTL2(m_ncoeffs > mode,
609  "calling argument mode is larger than total expansion order");
610 
611  for(i = 0; i < nquad1; ++i)
612  {
613  Vmath::Vcopy(nquad0,(NekDouble *)(base0.get() + mode0*nquad0),
614  1, &outarray[0]+i*nquad0,1);
615  }
616 
617  for(i = 0; i < nquad0; ++i)
618  {
619  Vmath::Vmul(nquad1,(NekDouble *)(base1.get() + mode1*nquad1),1,
620  &outarray[0]+i,nquad0,&outarray[0]+i,nquad0);
621  }
622  }
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:274
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:192

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

◆ 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.

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  }
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:617
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:376
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
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69

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

◆ 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.

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  GetTraceToElementMap(i,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  }
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:15
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:762
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:537
void GetTraceToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
Definition: StdExpansion.h:703
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:692
LibUtilities::NekManager< StdMatrixKey, DNekBlkMat, StdMatrixKey::opLess > m_stdStaticCondMatrixManager
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:265
std::shared_ptr< StdSegExp > StdSegExpSharedPtr
Definition: StdSegExp.h:46
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:372

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::StdRegions::StdExpansion::DetShapeType(), Blas::Dgemv(), Nektar::StdRegions::eMass, Nektar::StdRegions::StdExpansion::GetInteriorMap(), Nektar::StdRegions::StdExpansion::GetTraceToElementMap(), 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(), sign, Vmath::Vcopy(), and Vmath::Vsub().

◆ 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 1481 of file StdQuadExp.cpp.

1485  {
1486 
1488 
1489  if(inarray.get() == outarray.get())
1490  {
1491  Array<OneD,NekDouble> tmp(m_ncoeffs);
1492  Vmath::Vcopy(m_ncoeffs,inarray.get(),1,tmp.get(),1);
1493 
1494  Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, 1.0, mat->GetPtr().get(),
1495  m_ncoeffs, tmp.get(), 1, 0.0, outarray.get(), 1);
1496  }
1497  else
1498  {
1499  Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, 1.0, mat->GetPtr().get(),
1500  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
1501  }
1502  }
LibUtilities::NekManager< StdMatrixKey, DNekMat, StdMatrixKey::opLess > m_stdMatrixManager

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

◆ 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 1335 of file StdQuadExp.cpp.

1336  {
1337  int i;
1338  int order0 = GetBasisNumModes(0);
1339  int order1 = GetBasisNumModes(1);
1340  MatrixType mtype = mkey.GetMatrixType();
1341 
1342  DNekMatSharedPtr Mat;
1343 
1344  switch(mtype)
1345  {
1347  {
1348  int nq0 = m_base[0]->GetNumPoints();
1349  int nq1 = m_base[1]->GetNumPoints();
1350  int nq;
1351 
1352  // take definition from key
1353  if(mkey.ConstFactorExists(eFactorConst))
1354  {
1355  nq = (int) mkey.GetConstFactor(eFactorConst);
1356  }
1357  else
1358  {
1359  nq = max(nq0,nq1);
1360  }
1361 
1362  int neq = LibUtilities::StdQuadData::
1363  getNumberOfCoefficients(nq, nq);
1364  Array<OneD, Array<OneD, NekDouble> > coords(neq);
1365  Array<OneD, NekDouble> coll (2);
1366  Array<OneD, DNekMatSharedPtr> I (2);
1367  Array<OneD, NekDouble> tmp (nq0);
1368 
1370  AllocateSharedPtr(neq, nq0 * nq1);
1371  int cnt = 0;
1372 
1373  for(int i = 0; i < nq; ++i)
1374  {
1375  for(int j = 0; j < nq; ++j,++cnt)
1376  {
1377  coords[cnt] = Array<OneD, NekDouble>(2);
1378  coords[cnt][0] = -1.0 + 2*j/(NekDouble)(nq-1);
1379  coords[cnt][1] = -1.0 + 2*i/(NekDouble)(nq-1);
1380  }
1381  }
1382 
1383  for(int i = 0; i < neq; ++i)
1384  {
1385  LocCoordToLocCollapsed(coords[i],coll);
1386 
1387  I[0] = m_base[0]->GetI(coll);
1388  I[1] = m_base[1]->GetI(coll+1);
1389 
1390  // interpolate first coordinate direction
1391  for (int j = 0; j < nq1; ++j)
1392  {
1393  NekDouble fac = (I[1]->GetPtr())[j];
1394  Vmath::Smul(nq0,fac,I[0]->GetPtr(),1,tmp,1);
1395 
1396  Vmath::Vcopy(nq0, &tmp[0], 1,
1397  Mat->GetRawPtr()+j*nq0*neq+i,neq);
1398  }
1399 
1400  }
1401  break;
1402  }
1403  case eMass:
1404  {
1406  // For Fourier basis set the imaginary component of mean mode
1407  // to have a unit diagonal component in mass matrix
1409  {
1410  for(i = 0; i < order1; ++i)
1411  {
1412  (*Mat)(order0*i+1,i*order0+1) = 1.0;
1413  }
1414  }
1415 
1417  {
1418  for(i = 0; i < order0; ++i)
1419  {
1420  (*Mat)(order0+i ,order0+i) = 1.0;
1421  }
1422  }
1423  break;
1424  }
1425  case eFwdTrans:
1426  {
1429  StdMatrixKey iprodkey(eIProductWRTBase,DetShapeType(),*this);
1430  DNekMat &Iprod = *GetStdMatrix(iprodkey);
1431  StdMatrixKey imasskey(eInvMass,DetShapeType(),*this);
1432  DNekMat &Imass = *GetStdMatrix(imasskey);
1433 
1434  (*Mat) = Imass*Iprod;
1435  break;
1436  }
1437  case eGaussDG:
1438  {
1439  ConstFactorMap factors = mkey.GetConstFactors();
1440 
1441  int edge = (int)factors[StdRegions::eFactorGaussEdge];
1442  int dir = (edge + 1) % 2;
1443  int nCoeffs = m_base[dir]->GetNumModes();
1444 
1445  const LibUtilities::PointsKey BS_p(
1447  const LibUtilities::BasisKey BS_k(
1448  LibUtilities::eGauss_Lagrange, nCoeffs, BS_p);
1449 
1450  Array<OneD, NekDouble> coords(1, 0.0);
1451  coords[0] = (edge == 0 || edge == 3) ? -1.0 : 1.0;
1452 
1455  DNekMatSharedPtr m_Ix = basis->GetI(coords);
1456 
1458  1.0, nCoeffs);
1459  Vmath::Vcopy(nCoeffs, m_Ix->GetPtr(), 1, Mat->GetPtr(), 1);
1460  break;
1461  }
1462  default:
1463  {
1465  break;
1466  }
1467  }
1468 
1469  return Mat;
1470  }
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:158
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
Definition: StdExpansion.h:982
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:171
int getNumberOfCoefficients(int Na, int Nb)
Definition: ShapeType.hpp:137
BasisManagerT & BasisManager(void)
std::shared_ptr< Basis > BasisSharedPtr
@ eGaussGaussLegendre
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:48
@ eGauss_Lagrange
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:55
@ eFourier
Fourier Expansion .
Definition: BasisType.h:53
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:314
NekMatrix< NekDouble, StandardMatrixTag > DNekMat
Definition: NekTypeDefs.hpp:51
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:225

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

◆ v_GetBoundaryMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 798 of file StdQuadExp.cpp.

799  {
800  int i;
801  int cnt=0;
802  int nummodes0, nummodes1;
803  int value1 = 0, value2 = 0;
804  if(outarray.size()!=NumBndryCoeffs())
805  {
806  outarray = Array<OneD, unsigned int>(NumBndryCoeffs());
807  }
808 
809  nummodes0 = m_base[0]->GetNumModes();
810  nummodes1 = m_base[1]->GetNumModes();
811 
812  const LibUtilities::BasisType Btype0 = GetBasisType(0);
813  const LibUtilities::BasisType Btype1 = GetBasisType(1);
814 
815  switch(Btype1)
816  {
819  value1 = nummodes0;
820  break;
822  value1 = 2*nummodes0;
823  break;
824  default:
825  ASSERTL0(0,"Mapping array is not defined for this expansion");
826  break;
827  }
828 
829  for(i = 0; i < value1; i++)
830  {
831  outarray[i]=i;
832  }
833  cnt=value1;
834 
835  switch(Btype0)
836  {
839  value2 = value1+nummodes0-1;
840  break;
842  value2 = value1+1;
843  break;
844  default:
845  ASSERTL0(0,"Mapping array is not defined for this expansion");
846  break;
847  }
848 
849  for(i = 0; i < nummodes1-2; i++)
850  {
851  outarray[cnt++]=value1+i*nummodes0;
852  outarray[cnt++]=value2+i*nummodes0;
853  }
854 
855 
857  {
858  for(i = nummodes0*(nummodes1-1);i < GetNcoeffs(); i++)
859  {
860  outarray[cnt++] = i;
861  }
862  }
863  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:124
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:54
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48

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

◆ 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 744 of file StdQuadExp.cpp.

747  {
748  boost::ignore_unused(coords_2);
749  Array<OneD, const NekDouble> z0 = m_base[0]->GetZ();
750  Array<OneD, const NekDouble> z1 = m_base[1]->GetZ();
751  int nq0 = GetNumPoints(0);
752  int nq1 = GetNumPoints(1);
753  int i;
754 
755  for(i = 0; i < nq1; ++i)
756  {
757  Blas::Dcopy(nq0,z0.get(), 1,&coords_0[0] + i*nq0,1);
758  Vmath::Fill(nq0,z1[i],&coords_1[0] + i*nq0,1);
759  }
760  }
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:160
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45

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

◆ v_GetInteriorMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 865 of file StdQuadExp.cpp.

866  {
867  int i,j;
868  int cnt=0;
869  int nummodes0, nummodes1;
870  int startvalue = 0;
871  if(outarray.size()!=GetNcoeffs()-NumBndryCoeffs())
872  {
873  outarray = Array<OneD, unsigned int>(GetNcoeffs()-NumBndryCoeffs());
874  }
875 
876  nummodes0 = m_base[0]->GetNumModes();
877  nummodes1 = m_base[1]->GetNumModes();
878 
879  const LibUtilities::BasisType Btype0 = GetBasisType(0);
880  const LibUtilities::BasisType Btype1 = GetBasisType(1);
881 
882  switch(Btype1)
883  {
885  startvalue = nummodes0;
886  break;
888  startvalue = 2*nummodes0;
889  break;
890  default:
891  ASSERTL0(0,"Mapping array is not defined for this expansion");
892  break;
893  }
894 
895  switch(Btype0)
896  {
898  startvalue++;
899  break;
901  startvalue+=2;
902  break;
903  default:
904  ASSERTL0(0,"Mapping array is not defined for this expansion");
905  break;
906  }
907 
908  for(i = 0; i < nummodes1-2; i++)
909  {
910  for(j = 0; j < nummodes0-2; j++)
911  {
912  outarray[cnt++]=startvalue+j;
913  }
914  startvalue+=nummodes0;
915  }
916  }

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

◆ v_GetNtraces()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 633 of file StdQuadExp.cpp.

634  {
635  return 4;
636  }

◆ v_GetNverts()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 628 of file StdQuadExp.cpp.

629  {
630  return 4;
631  }

◆ v_GetSimplexEquiSpacedConnectivity()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1772 of file StdQuadExp.cpp.

1775  {
1776  boost::ignore_unused(standard);
1777 
1778  int np1 = m_base[0]->GetNumPoints();
1779  int np2 = m_base[1]->GetNumPoints();
1780  int np = max(np1,np2);
1781 
1782  conn = Array<OneD, int>(6*(np-1)*(np-1));
1783 
1784  int row = 0;
1785  int rowp1 = 0;
1786  int cnt = 0;
1787  for(int i = 0; i < np-1; ++i)
1788  {
1789  rowp1 += np;
1790  for(int j = 0; j < np-1; ++j)
1791  {
1792  conn[cnt++] = row +j;
1793  conn[cnt++] = row +j+1;
1794  conn[cnt++] = rowp1 +j;
1795 
1796  conn[cnt++] = rowp1 +j+1;
1797  conn[cnt++] = rowp1 +j;
1798  conn[cnt++] = row +j+1;
1799  }
1800  row += np;
1801  }
1802  }

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

◆ v_GetTraceBasisKey()

const LibUtilities::BasisKey Nektar::StdRegions::StdQuadExp::v_GetTraceBasisKey ( const int  i,
const int  j 
) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 666 of file StdQuadExp.cpp.

668  {
669  boost::ignore_unused(j);
670  ASSERTL2((i >= 0)&&(i <= 3),"edge id is out of range");
671 
672  if((i == 0)||(i == 2))
673  {
674  return GetBasis(0)->GetBasisKey();
675  }
676  else
677  {
678  return GetBasis(1)->GetBasisKey();
679  }
680 
681  }
const LibUtilities::BasisSharedPtr & GetBasis(int dir) const
This function gets the shared point to basis in the dir direction.
Definition: StdExpansion.h:111

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

◆ v_GetTraceInteriorToElementMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1204 of file StdQuadExp.cpp.

1209  {
1210  int i;
1211  const int nummodes0 = m_base[0]->GetNumModes();
1212  const int nummodes1 = m_base[1]->GetNumModes();
1213  const int nEdgeIntCoeffs = GetTraceNcoeffs(eid)-2;
1214  const LibUtilities::BasisType bType = GetBasisType(eid%2);
1215 
1216  if(maparray.size() != nEdgeIntCoeffs)
1217  {
1218  maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
1219  }
1220 
1221  if(signarray.size() != nEdgeIntCoeffs)
1222  {
1223  signarray = Array<OneD, int>(nEdgeIntCoeffs,1);
1224  }
1225  else
1226  {
1227  fill( signarray.get() , signarray.get()+nEdgeIntCoeffs, 1 );
1228  }
1229 
1230  if(bType == LibUtilities::eModified_A)
1231  {
1232  switch(eid)
1233  {
1234  case 0:
1235  {
1236  for(i = 0; i < nEdgeIntCoeffs; i++)
1237  {
1238  maparray[i] = i+2;
1239  }
1240  }
1241  break;
1242  case 1:
1243  {
1244  for(i = 0; i < nEdgeIntCoeffs; i++)
1245  {
1246  maparray[i] = (i+2)*nummodes0 + 1;
1247  }
1248  }
1249  break;
1250  case 2:
1251  {
1252  for(i = 0; i < nEdgeIntCoeffs; i++)
1253  {
1254  maparray[i] = nummodes0+i+2;
1255  }
1256  }
1257  break;
1258  case 3:
1259  {
1260  for(i = 0; i < nEdgeIntCoeffs; i++)
1261  {
1262  maparray[i] = (i+2)*nummodes0;
1263  }
1264  }
1265  break;
1266  default:
1267  ASSERTL0(false,"eid must be between 0 and 3");
1268  break;
1269  }
1270 
1271  if(edgeOrient==eBackwards)
1272  {
1273  for(i = 1; i < nEdgeIntCoeffs; i+=2)
1274  {
1275  signarray[i] = -1;
1276  }
1277  }
1278  }
1279  else if(bType == LibUtilities::eGLL_Lagrange)
1280  {
1281  switch(eid)
1282  {
1283  case 0:
1284  {
1285  for(i = 0; i < nEdgeIntCoeffs; i++)
1286  {
1287  maparray[i] = i+1;
1288  }
1289  }
1290  break;
1291  case 1:
1292  {
1293  for(i = 0; i < nEdgeIntCoeffs; i++)
1294  {
1295  maparray[i] = (i+2)*nummodes0 - 1;
1296  }
1297  }
1298  break;
1299  case 2:
1300  {
1301  for(i = 0; i < nEdgeIntCoeffs; i++)
1302  {
1303  maparray[i] = nummodes0*(nummodes1-1) + i + 1;
1304  }
1305  }
1306  break;
1307  case 3:
1308  {
1309  for(i = 0; i < nEdgeIntCoeffs; i++)
1310  {
1311  maparray[i] = nummodes0 * (i+1);
1312  }
1313  }
1314  break;
1315  default:
1316  ASSERTL0(false,"eid must be between 0 and 3");
1317  break;
1318  }
1319  if(edgeOrient == eBackwards)
1320  {
1321  reverse( maparray.get() , maparray.get()+nEdgeIntCoeffs );
1322  }
1323  }
1324  else
1325  {
1326  ASSERTL0(false,"Mapping not defined for this type of basis");
1327  }
1328 
1329  }
int GetTraceNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th trace.
Definition: StdExpansion.h:265

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

◆ v_GetTraceNcoeffs()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 638 of file StdQuadExp.cpp.

639  {
640  ASSERTL2((i >= 0)&&(i <= 3),"edge id is out of range");
641 
642  if((i == 0)||(i == 2))
643  {
644  return GetBasisNumModes(0);
645  }
646  else
647  {
648  return GetBasisNumModes(1);
649  }
650  }

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

◆ v_GetTraceNumPoints()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 652 of file StdQuadExp.cpp.

653  {
654  ASSERTL2((i >= 0)&&(i <= 3),"edge id is out of range");
655 
656  if((i == 0)||(i == 2))
657  {
658  return GetNumPoints(0);
659  }
660  else
661  {
662  return GetNumPoints(1);
663  }
664  }

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

◆ v_GetTraceToElementMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1026 of file StdQuadExp.cpp.

1032  {
1033  boost::ignore_unused(Q);
1034  int i;
1035  int numModes=0;
1036  int order0 = m_base[0]->GetNumModes();
1037  int order1 = m_base[1]->GetNumModes();
1038 
1039  switch (eid)
1040  {
1041  case 0:
1042  case 2:
1043  numModes = order0;
1044  break;
1045  case 1:
1046  case 3:
1047  numModes = order1;
1048  break;
1049  default:
1050  ASSERTL0(false,"eid must be between 0 and 3");
1051  }
1052 
1053  bool checkForZeroedModes = false;
1054  if (P == -1)
1055  {
1056  P = numModes;
1057  }
1058  else if(P != numModes)
1059  {
1060  checkForZeroedModes = true;
1061  }
1062 
1063 
1064  if (maparray.size() != P)
1065  {
1066  maparray = Array<OneD, unsigned int>(P);
1067  }
1068 
1069  if(signarray.size() != P)
1070  {
1071  signarray = Array<OneD, int>(P, 1);
1072  }
1073  else
1074  {
1075  fill(signarray.get(), signarray.get()+P, 1);
1076  }
1077 
1078  const LibUtilities::BasisType bType = GetBasisType(eid%2);
1079 
1080  if (bType == LibUtilities::eModified_A)
1081  {
1082  switch (eid)
1083  {
1084  case 0:
1085  {
1086  for (i = 0; i < P; i++)
1087  {
1088  maparray[i] = i;
1089  }
1090  }
1091  break;
1092  case 1:
1093  {
1094  for (i = 0; i < P; i++)
1095  {
1096  maparray[i] = i*order0 + 1;
1097  }
1098  }
1099  break;
1100  case 2:
1101  {
1102  for (i = 0; i < P; i++)
1103  {
1104  maparray[i] = order0+i;
1105  }
1106  }
1107  break;
1108  case 3:
1109  {
1110  for (i = 0; i < P; i++)
1111  {
1112  maparray[i] = i*order0;
1113  }
1114  }
1115  break;
1116  default:
1117  ASSERTL0(false, "eid must be between 0 and 3");
1118  break;
1119  }
1120 
1121  if (edgeOrient == eBackwards)
1122  {
1123  swap(maparray[0], maparray[1]);
1124 
1125  for(i = 3; i < P; i+=2)
1126  {
1127  signarray[i] = -1;
1128  }
1129  }
1130  }
1131  else if(bType == LibUtilities::eGLL_Lagrange ||
1133  {
1134  switch (eid)
1135  {
1136  case 0:
1137  {
1138  for (i = 0; i < P; i++)
1139  {
1140  maparray[i] = i;
1141  }
1142  }
1143  break;
1144  case 1:
1145  {
1146  for (i = 0; i < P; i++)
1147  {
1148  maparray[i] = (i+1)*order0 - 1;
1149  }
1150  }
1151  break;
1152  case 2:
1153  {
1154  for (i = 0; i < P; i++)
1155  {
1156  maparray[i] = order0*(order1-1) + i;
1157  }
1158  }
1159  break;
1160  case 3:
1161  {
1162  for (i = 0; i < P; i++)
1163  {
1164  maparray[i] = order0*i;
1165  }
1166  }
1167  break;
1168  default:
1169  ASSERTL0(false, "eid must be between 0 and 3");
1170  break;
1171  }
1172 
1173  if (edgeOrient == eBackwards)
1174  {
1175  reverse(maparray.get(), maparray.get()+P);
1176  }
1177  }
1178  else
1179  {
1180  ASSERTL0(false, "Mapping not defined for this type of basis");
1181  }
1182 
1183  if (checkForZeroedModes)
1184  {
1185  if (bType == LibUtilities::eModified_A)
1186  {
1187  // Zero signmap and set maparray to zero if
1188  // elemental modes are not as large as face modes
1189  for (int j = numModes; j < P; j++)
1190  {
1191  signarray[j] = 0.0;
1192  maparray[j] = maparray[0];
1193  }
1194  }
1195  else
1196  {
1197  ASSERTL0(false, "Different trace space edge dimension "
1198  "and element edge dimension not possible "
1199  "for GLL-Lagrange bases");
1200  }
1201  }
1202  }

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

◆ v_GetVertexMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 918 of file StdQuadExp.cpp.

919  {
920  int localDOF = 0;
921 
922  if(useCoeffPacking == true)
923  {
924  switch(localVertexId)
925  {
926  case 0:
927  {
928  localDOF = 0;
929  }
930  break;
931  case 1:
932  {
934  {
935  localDOF = m_base[0]->GetNumModes()-1;
936  }
937  else
938  {
939  localDOF = 1;
940  }
941  }
942  break;
943  case 2:
944  {
946  {
947  localDOF = m_base[0]->GetNumModes() * (m_base[1]->GetNumModes()-1);
948  }
949  else
950  {
951  localDOF = m_base[0]->GetNumModes();
952  }
953  }
954  break;
955  case 3:
956  {
958  {
959  localDOF = m_base[0]->GetNumModes()*m_base[1]->GetNumModes()-1;
960  }
961  else
962  {
963  localDOF = m_base[0]->GetNumModes()+1;
964  }
965  }
966  break;
967  default:
968  ASSERTL0(false,"eid must be between 0 and 3");
969  break;
970  }
971 
972  }
973  else
974  {
975  switch(localVertexId)
976  {
977  case 0:
978  {
979  localDOF = 0;
980  }
981  break;
982  case 1:
983  {
985  {
986  localDOF = m_base[0]->GetNumModes()-1;
987  }
988  else
989  {
990  localDOF = 1;
991  }
992  }
993  break;
994  case 2:
995  {
997  {
998  localDOF = m_base[0]->GetNumModes()*m_base[1]->GetNumModes()-1;
999  }
1000  else
1001  {
1002  localDOF = m_base[0]->GetNumModes()+1;
1003  }
1004  }
1005  break;
1006  case 3:
1007  {
1009  {
1010  localDOF = m_base[0]->GetNumModes() * (m_base[1]->GetNumModes()-1);
1011  }
1012  else
1013  {
1014  localDOF = m_base[0]->GetNumModes();
1015  }
1016  }
1017  break;
1018  default:
1019  ASSERTL0(false,"eid must be between 0 and 3");
1020  break;
1021  }
1022  }
1023  return localDOF;
1024  }

References ASSERTL0, Nektar::LibUtilities::eGLL_Lagrange, Nektar::StdRegions::StdExpansion::GetBasisType(), and Nektar::StdRegions::StdExpansion::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 1739 of file StdQuadExp.cpp.

1742  {
1743  StdQuadExp::v_HelmholtzMatrixOp_MatFree(inarray,outarray,mkey);
1744  }
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)

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

◆ 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.

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)

References Nektar::StdRegions::StdExpansion2D::Integral(), and Nektar::StdRegions::StdExpansion::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_p \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.

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:733
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
Definition: StdQuadExp.cpp:399

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

Referenced by v_FwdTrans().

◆ 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.

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  }
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:134

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

◆ 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.

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

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

Referenced by v_IProductWRTBase().

◆ 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.

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.size()>=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 double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = .
Definition: Blas.hpp:197

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

◆ 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.

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)

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

◆ 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.

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  }

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.

◆ 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.

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  }

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

◆ v_IsBoundaryInteriorExpansion()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 727 of file StdQuadExp.cpp.

728  {
729  bool returnval = false;
730 
733  {
736  {
737  returnval = true;
738  }
739  }
740 
741  return returnval;
742  }

References Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), and Nektar::StdRegions::StdExpansion::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 1715 of file StdQuadExp.cpp.

1719  {
1720  StdQuadExp::v_LaplacianMatrixOp_MatFree(inarray,outarray,mkey);
1721  }
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)

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

◆ 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 1723 of file StdQuadExp.cpp.

1727  {
1728  StdExpansion::LaplacianMatrixOp_MatFree(k1,k2,inarray,outarray,mkey);
1729  }
void LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)

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

◆ v_LocCollapsedToLocCoord()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 580 of file StdQuadExp.cpp.

582  {
583  xi[0] = eta[0];
584  xi[1] = eta[1];
585  }

◆ 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 1707 of file StdQuadExp.cpp.

1711  {
1712  StdExpansion::MassMatrixOp_MatFree(inarray,outarray,mkey);
1713  }
void MassMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)

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

◆ 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 1747 of file StdQuadExp.cpp.

1750  {
1751  int i;
1752  int nquad0 = m_base[0]->GetNumPoints();
1753  int nquad1 = m_base[1]->GetNumPoints();
1754 
1755  const Array<OneD, const NekDouble>& w0 = m_base[0]->GetW();
1756  const Array<OneD, const NekDouble>& w1 = m_base[1]->GetW();
1757 
1758  // multiply by integration constants
1759  for(i = 0; i < nquad1; ++i)
1760  {
1761  Vmath::Vmul(nquad0, inarray.get()+i*nquad0,1,
1762  w0.get(),1,outarray.get()+i*nquad0,1);
1763  }
1764 
1765  for(i = 0; i < nquad0; ++i)
1766  {
1767  Vmath::Vmul(nquad1,outarray.get()+i,nquad0,w1.get(),1,
1768  outarray.get()+i,nquad0);
1769  }
1770  }

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

◆ v_NumBndryCoeffs()

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

◆ v_NumDGBndryCoeffs()

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

◆ 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.

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.

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

Referenced by v_StdPhysDeriv().

◆ 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.

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

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

◆ v_PhysEvaluateBasis()

NekDouble Nektar::StdRegions::StdQuadExp::v_PhysEvaluateBasis ( const Array< OneD, const NekDouble > &  coords,
int  mode 
)
protectedvirtual

This function evaluates the basis function mode mode at a point coords of the domain.

This function uses barycentric interpolation with the tensor product separation of the basis function to improve performance.

Parameters
coordThe coordinate inside the standard region.
modeThe mode number to be evaluated.
Returns
The value of the basis function mode at coords.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 774 of file StdQuadExp.cpp.

777  {
778  ASSERTL2(coords[0] > -1 - NekConstants::kNekZeroTol,
779  "coord[0] < -1");
780  ASSERTL2(coords[0] < 1 + NekConstants::kNekZeroTol,
781  "coord[0] > 1");
782  ASSERTL2(coords[1] > -1 - NekConstants::kNekZeroTol,
783  "coord[1] < -1");
784  ASSERTL2(coords[1] < 1 + NekConstants::kNekZeroTol,
785  "coord[1] > 1");
786 
787  const int nm0 = m_base[0]->GetNumModes();
788  const int nm1 = m_base[1]->GetNumModes();
789 
790  return StdExpansion::BaryEvaluateBasis<0>(coords[0], mode % nm1) *
791  StdExpansion::BaryEvaluateBasis<1>(coords[1], mode / nm0);
792  }
static const NekDouble kNekZeroTol

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

◆ 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 1657 of file StdQuadExp.cpp.

1661  {
1662  int n_coeffs = inarray.size();
1663 
1664 
1665  Array<OneD, NekDouble> coeff(n_coeffs);
1666  Array<OneD, NekDouble> coeff_tmp(n_coeffs,0.0);
1667  Array<OneD, NekDouble> tmp;
1668  Array<OneD, NekDouble> tmp2;
1669 
1670  int nmodes0 = m_base[0]->GetNumModes();
1671  int nmodes1 = m_base[1]->GetNumModes();
1672  int numMax = nmodes0;
1673 
1674  Vmath::Vcopy(n_coeffs,inarray,1,coeff_tmp,1);
1675 
1676  const LibUtilities::PointsKey Pkey0(
1678  const LibUtilities::PointsKey Pkey1(
1680 
1681  LibUtilities::BasisKey b0(m_base[0]->GetBasisType(),nmodes0,Pkey0);
1682  LibUtilities::BasisKey b1(m_base[1]->GetBasisType(),nmodes1,Pkey1);
1683 
1684  LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A,nmodes0,Pkey0);
1685  LibUtilities::BasisKey bortho1(LibUtilities::eOrtho_A,nmodes1,Pkey1);
1686 
1688  b0, b1, coeff_tmp, bortho0, bortho1, coeff);
1689 
1690  Vmath::Zero(n_coeffs,coeff_tmp,1);
1691 
1692  int cnt = 0;
1693  for (int i = 0; i < numMin+1; ++i)
1694  {
1695  Vmath::Vcopy(numMin,
1696  tmp = coeff+cnt,1,
1697  tmp2 = coeff_tmp+cnt,1);
1698 
1699  cnt = i*numMax;
1700  }
1701 
1703  bortho0, bortho1, coeff_tmp, b0, b1, outarray);
1704  }
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
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:436

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

◆ 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.

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

References v_PhysDeriv().

◆ 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.

145  {
146  //PhysTensorDeriv(inarray, outarray);
147  StdQuadExp::v_PhysDeriv(dir,inarray,outarray);
148 
149  }

References v_PhysDeriv().

◆ 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 1505 of file StdQuadExp.cpp.

1507  {
1508  // Generate an orthonogal expansion
1509  int qa = m_base[0]->GetNumPoints();
1510  int qb = m_base[1]->GetNumPoints();
1511  int nmodes_a = m_base[0]->GetNumModes();
1512  int nmodes_b = m_base[1]->GetNumModes();
1513  int nmodes = min(nmodes_a,nmodes_b);
1514  // Declare orthogonal basis.
1515  LibUtilities::PointsKey pa(qa,m_base[0]->GetPointsType());
1516  LibUtilities::PointsKey pb(qb,m_base[1]->GetPointsType());
1517 
1518  LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A,nmodes_a,pa);
1519  LibUtilities::BasisKey Bb(LibUtilities::eOrtho_A,nmodes_b,pb);
1520  StdQuadExp OrthoExp(Ba,Bb);
1521 
1522  //SVV parameters loaded from the .xml case file
1523  Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
1524 
1525  // project onto modal space.
1526  OrthoExp.FwdTrans(array,orthocoeffs);
1527 
1528  if(mkey.ConstFactorExists(eFactorSVVPowerKerDiffCoeff)) // Rodrigo's power kernel
1529  {
1530  NekDouble cutoff = mkey.GetConstFactor(eFactorSVVCutoffRatio);
1531  NekDouble SvvDiffCoeff =
1532  mkey.GetConstFactor(eFactorSVVPowerKerDiffCoeff)*
1533  mkey.GetConstFactor(eFactorSVVDiffCoeff);
1534 
1535  for(int j = 0; j < nmodes_a; ++j)
1536  {
1537  for(int k = 0; k < nmodes_b; ++k)
1538  {
1539  // linear space but makes high modes very negative
1540  NekDouble fac = std::max(
1541  pow((1.0*j)/(nmodes_a-1),cutoff*nmodes_a),
1542  pow((1.0*k)/(nmodes_b-1),cutoff*nmodes_b));
1543 
1544  orthocoeffs[j*nmodes_b+k] *= SvvDiffCoeff * fac;
1545  }
1546  }
1547  }
1548  else if(mkey.ConstFactorExists(eFactorSVVDGKerDiffCoeff)) // Rodrigo/mansoor's DG kernel
1549  {
1550  NekDouble SvvDiffCoeff =
1551  mkey.GetConstFactor(eFactorSVVDGKerDiffCoeff)*
1552  mkey.GetConstFactor(eFactorSVVDiffCoeff);
1553  int max_ab = max(nmodes_a-kSVVDGFiltermodesmin,
1554  nmodes_b-kSVVDGFiltermodesmin);
1555  max_ab = max(max_ab,0);
1556  max_ab = min(max_ab,kSVVDGFiltermodesmax-kSVVDGFiltermodesmin);
1557 
1558  for(int j = 0; j < nmodes_a; ++j)
1559  {
1560  for(int k = 0; k < nmodes_b; ++k)
1561  {
1562  int maxjk = max(j,k);
1563  maxjk = min(maxjk,kSVVDGFiltermodesmax-1);
1564 
1565  orthocoeffs[j*nmodes_b+k] *= SvvDiffCoeff *
1566  kSVVDGFilter[max_ab][maxjk];
1567  }
1568  }
1569  }
1570  else
1571  {
1572  NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDiffCoeff);
1573  //Exponential Kernel implementation
1574  int cutoff = (int) (mkey.GetConstFactor(eFactorSVVCutoffRatio)*
1575  min(nmodes_a,nmodes_b));
1576 
1577  //counters for scanning through orthocoeffs array
1578  int cnt = 0;
1579 
1580  //------"New" Version August 22nd '13--------------------
1581  for(int j = 0; j < nmodes_a; ++j)
1582  {
1583  for(int k = 0; k < nmodes_b; ++k)
1584  {
1585  if(j + k >= cutoff)//to filter out only the "high-modes"
1586  {
1587  orthocoeffs[j*nmodes_b+k] *=
1588  (SvvDiffCoeff*exp(-(j+k-nmodes)*(j+k-nmodes)/
1589  ((NekDouble)((j+k-cutoff+1)*
1590  (j+k-cutoff+1)))));
1591  }
1592  else
1593  {
1594  orthocoeffs[j*nmodes_b+k] *= 0.0;
1595  }
1596  cnt++;
1597  }
1598  }
1599  }
1600 
1601  // backward transform to physical space
1602  OrthoExp.BwdTrans(orthocoeffs,array);
1603  }
const int kSVVDGFiltermodesmin
Definition: StdRegions.hpp:388
const int kSVVDGFiltermodesmax
Definition: StdRegions.hpp:389
const NekDouble kSVVDGFilter[9][11]
Definition: StdRegions.hpp:391

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.

◆ 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 1731 of file StdQuadExp.cpp.

1735  {
1736  StdExpansion::WeakDerivMatrixOp_MatFree(i,inarray,outarray,mkey);
1737  }
void WeakDerivMatrixOp_MatFree(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)

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