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

#include <StdQuadExp.h>

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

Public Member Functions

 StdQuadExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb)
 Constructor using BasisKey class for quadrature points and order definition. More...
 
 StdQuadExp (const StdQuadExp &T)=default
 
 ~StdQuadExp () override=default
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion2D
 StdExpansion2D (int numcoeffs, const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb)
 
 StdExpansion2D ()=default
 
 StdExpansion2D (const StdExpansion2D &T)=default
 
 ~StdExpansion2D () override=default
 
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)
 
NekDouble BaryTensorDeriv (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
 
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 () const
 
std::shared_ptr< StdExpansionGetLinStdExp (void) const
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion () const
 
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 FwdTransBndConstrained (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 GetTraceCoeffMap (const unsigned int traceid, Array< OneD, unsigned int > &maparray)
 
void GetElmtTraceToTraceMap (const unsigned 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 LinearAdvectionMatrixOp (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, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
 This function evaluates the first derivative of the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs, std::array< NekDouble, 6 > &secondOrderDerivs)
 
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 int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
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) override
 Integrates the specified function over the domain. More...
 
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) override
 Calculate the derivative of the physical points. More...
 
void v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Calculate the derivative of the physical points in a given direction. More...
 
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) override
 
void v_StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
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) override
 
void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Transform a given function from physical quadrature space to coefficient space. More...
 
void v_FwdTransBndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Calculate the inner product of inarray with respect to the basis B=base0*base1 and put into outarray. More...
 
void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
 
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) override
 
void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta) override
 
void v_LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi) override
 
void v_FillMode (const int mode, Array< OneD, NekDouble > &array) override
 Fill outarray with mode mode of expansion. More...
 
int v_GetNverts () const final
 
int v_GetNtraces () const final
 
int v_GetTraceNcoeffs (const int i) const final
 
int v_GetTraceIntNcoeffs (const int i) const final
 
int v_GetTraceNumPoints (const int i) const final
 
int v_NumBndryCoeffs () const final
 
int v_NumDGBndryCoeffs () const final
 
int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset) override
 
const LibUtilities::BasisKey v_GetTraceBasisKey (const int i, const int j) const final
 
LibUtilities::ShapeType v_DetShapeType () const final
 
bool v_IsBoundaryInteriorExpansion () const override
 
void v_GetCoords (Array< OneD, NekDouble > &coords_0, Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2) override
 
NekDouble v_PhysEvaluateBasis (const Array< OneD, const NekDouble > &coords, int mode) override
 This function evaluates the basis function mode mode at a point coords of the domain. More...
 
NekDouble v_PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
 
void v_GetBoundaryMap (Array< OneD, unsigned int > &outarray) override
 
void v_GetInteriorMap (Array< OneD, unsigned int > &outarray) override
 
int v_GetVertexMap (int localVertexId, bool useCoeffPacking=false) override
 
void v_GetTraceCoeffMap (const unsigned int traceid, Array< OneD, unsigned int > &maparray) override
 Get the map of the coefficient location to teh local trace coefficients. More...
 
void v_GetTraceInteriorToElementMap (const int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation edgeOrient=eForwards) override
 
DNekMatSharedPtr v_GenMatrix (const StdMatrixKey &mkey) override
 
DNekMatSharedPtr v_CreateStdMatrix (const StdMatrixKey &mkey) override
 
void v_MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
 
void v_LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
 
void v_LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
 
void v_WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
 
void v_HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
 
void v_SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey) override
 
void v_ExponentialFilter (Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff) override
 
void v_ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true) override
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion2D
NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals) override
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble v_PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals) override
 
NekDouble v_PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
 
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)=0
 
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)=0
 
void v_LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
void v_HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
void v_GetTraceCoeffMap (const unsigned int traceid, Array< OneD, unsigned int > &maparray) override
 
void v_GetElmtTraceToTraceMap (const unsigned int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient, int P, int Q) override
 Determine the mapping to re-orientate the coefficients along the element trace (assumed to align with the standard element) into the orientation of the local trace given by edgeOrient. More...
 
void v_GetTraceToElementMap (const int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient=eForwards, int P=-1, int Q=-1) override
 
void v_GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat) override
 
- 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_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 LinearAdvectionMatrixOp_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)
 
virtual void v_GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 
virtual void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals, NekDouble &deriv, NekDouble &deriv2)
 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)
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals)
 Helper function to pass an unused value by reference into BaryEvaluate. More...
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals, NekDouble &deriv)
 

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 44 of file StdQuadExp.h.

Constructor & Destructor Documentation

◆ StdQuadExp() [1/2]

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

49 : StdExpansion(Ba.GetNumModes() * Bb.GetNumModes(), 2, Ba, Bb),
50 StdExpansion2D(Ba.GetNumModes() * Bb.GetNumModes(), Ba, Bb)
51{
52}
StdExpansion()
Default Constructor.

◆ StdQuadExp() [2/2]

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

◆ ~StdQuadExp()

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

Member Function Documentation

◆ v_BwdTrans()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 127 of file StdQuadExp.cpp.

129{
130 if (m_base[0]->Collocation() && m_base[1]->Collocation())
131 {
133 inarray, 1, outarray, 1);
134 }
135 else
136 {
137 StdQuadExp::v_BwdTrans_SumFac(inarray, outarray);
138 }
139}
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:218
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdQuadExp.cpp:141
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 141 of file StdQuadExp.cpp.

143{
144 Array<OneD, NekDouble> wsp(m_base[0]->GetNumPoints() *
145 m_base[1]->GetNumModes());
146
147 BwdTrans_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetBdata(), inarray,
148 outarray, wsp, true, true);
149}
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 
)
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion2D.

Definition at line 161 of file StdQuadExp.cpp.

167{
168 int nquad0 = m_base[0]->GetNumPoints();
169 int nquad1 = m_base[1]->GetNumPoints();
170 int nmodes0 = m_base[0]->GetNumModes();
171 int nmodes1 = m_base[1]->GetNumModes();
172
173 bool colldir0 = doCheckCollDir0 ? (m_base[0]->Collocation()) : false;
174 bool colldir1 = doCheckCollDir1 ? (m_base[1]->Collocation()) : false;
175
176 if (colldir0 && colldir1)
177 {
178 Vmath::Vcopy(m_ncoeffs, inarray.get(), 1, outarray.get(), 1);
179 }
180 else if (colldir0)
181 {
182 Blas::Dgemm('N', 'T', nquad0, nquad1, nmodes1, 1.0, &inarray[0], nquad0,
183 base1.get(), nquad1, 0.0, &outarray[0], nquad0);
184 }
185 else if (colldir1)
186 {
187 Blas::Dgemm('N', 'N', nquad0, nmodes1, nmodes0, 1.0, base0.get(),
188 nquad0, &inarray[0], nmodes0, 0.0, &outarray[0], nquad0);
189 }
190 else
191 {
192 ASSERTL1(wsp.size() >= nquad0 * nmodes1,
193 "Workspace size is not sufficient");
194
195 // Those two calls correpsond to the operation
196 // out = B0*in*Transpose(B1);
197 Blas::Dgemm('N', 'N', nquad0, nmodes1, nmodes0, 1.0, base0.get(),
198 nquad0, &inarray[0], nmodes0, 0.0, &wsp[0], nquad0);
199 Blas::Dgemm('N', 'T', nquad0, nquad1, nmodes1, 1.0, &wsp[0], nquad0,
200 base1.get(), nquad1, 0.0, &outarray[0], nquad0);
201 }
202}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
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:383

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 654 of file StdQuadExp.cpp.

656{
657 int nmodes = nummodes[modes_offset] * nummodes[modes_offset + 1];
658 modes_offset += 2;
659
660 return nmodes;
661}

◆ v_CreateStdMatrix()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 1329 of file StdQuadExp.cpp.

1330{
1331 return GenMatrix(mkey);
1332}
DNekMatSharedPtr GenMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:844

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

◆ v_DetShapeType()

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

◆ v_ExponentialFilter()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1436 of file StdQuadExp.cpp.

1440{
1441 // Generate an orthogonal expansion
1442 int qa = m_base[0]->GetNumPoints();
1443 int qb = m_base[1]->GetNumPoints();
1444 int nmodesA = m_base[0]->GetNumModes();
1445 int nmodesB = m_base[1]->GetNumModes();
1446 int P = nmodesA - 1;
1447 int Q = nmodesB - 1;
1448
1449 // Declare orthogonal basis.
1450 LibUtilities::PointsKey pa(qa, m_base[0]->GetPointsType());
1451 LibUtilities::PointsKey pb(qb, m_base[1]->GetPointsType());
1452
1453 LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A, nmodesA, pa);
1454 LibUtilities::BasisKey Bb(LibUtilities::eOrtho_A, nmodesB, pb);
1455 StdQuadExp OrthoExp(Ba, Bb);
1456
1457 // Cutoff
1458 int Pcut = cutoff * P;
1459 int Qcut = cutoff * Q;
1460
1461 // Project onto orthogonal space.
1462 Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
1463 OrthoExp.FwdTrans(array, orthocoeffs);
1464
1465 //
1466 NekDouble fac, fac1, fac2;
1467 for (int i = 0; i < nmodesA; ++i)
1468 {
1469 for (int j = 0; j < nmodesB; ++j)
1470 {
1471 // to filter out only the "high-modes"
1472 if (i > Pcut || j > Qcut)
1473 {
1474 fac1 = (NekDouble)(i - Pcut) / ((NekDouble)(P - Pcut));
1475 fac2 = (NekDouble)(j - Qcut) / ((NekDouble)(Q - Qcut));
1476 fac = max(fac1, fac2);
1477 fac = pow(fac, exponent);
1478 orthocoeffs[i * nmodesB + j] *= exp(-alpha * fac);
1479 }
1480 }
1481 }
1482
1483 // backward transform to physical space
1484 OrthoExp.BwdTrans(orthocoeffs, array);
1485}
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:205
StdQuadExp(const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb)
Constructor using BasisKey class for quadrature points and order definition.
Definition: StdQuadExp.cpp:47
@ P
Monomial polynomials .
Definition: BasisType.h:62
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:42
double NekDouble

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 Nektar::LibUtilities::P.

◆ v_FillMode()

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

Fill outarray with mode mode of expansion.

Note for quadrilateral expansions _base[0] (i.e. p) modes run fastest

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 521 of file StdQuadExp.cpp.

522{
523 int i;
524 int nquad0 = m_base[0]->GetNumPoints();
525 int nquad1 = m_base[1]->GetNumPoints();
526 Array<OneD, const NekDouble> base0 = m_base[0]->GetBdata();
527 Array<OneD, const NekDouble> base1 = m_base[1]->GetBdata();
528 int btmp0 = m_base[0]->GetNumModes();
529 int mode0 = mode % btmp0;
530 int mode1 = mode / btmp0;
531
532 ASSERTL2(mode1 == (int)floor((1.0 * mode) / btmp0),
533 "Integer Truncation not Equiv to Floor");
534
535 ASSERTL2(m_ncoeffs > mode,
536 "calling argument mode is larger than total expansion order");
537
538 for (i = 0; i < nquad1; ++i)
539 {
540 Vmath::Vcopy(nquad0, (NekDouble *)(base0.get() + mode0 * nquad0), 1,
541 &outarray[0] + i * nquad0, 1);
542 }
543
544 for (i = 0; i < nquad0; ++i)
545 {
546 Vmath::Vmul(nquad1, (NekDouble *)(base1.get() + mode1 * nquad1), 1,
547 &outarray[0] + i, nquad0, &outarray[0] + i, nquad0);
548 }
549}
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:265
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.hpp:72

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

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

206{
207 if ((m_base[0]->Collocation()) && (m_base[1]->Collocation()))
208 {
209 Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
210 }
211 else
212 {
213 StdQuadExp::v_IProductWRTBase(inarray, outarray);
214
215 // get Mass matrix inverse
216 StdMatrixKey masskey(eInvMass, DetShapeType(), *this);
217 DNekMatSharedPtr matsys = GetStdMatrix(masskey);
218
219 // copy inarray in case inarray == outarray
220 NekVector<NekDouble> in(m_ncoeffs, outarray, eCopy);
221 NekVector<NekDouble> out(m_ncoeffs, outarray, eWrapper);
222
223 out = (*matsys) * in;
224 }
225}
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:603
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:367
void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Calculate the inner product of inarray with respect to the basis B=base0*base1 and put into outarray.
Definition: StdQuadExp.cpp:354
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 227 of file StdQuadExp.cpp.

230{
231 if ((m_base[0]->Collocation()) && (m_base[1]->Collocation()))
232 {
233 Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
234 }
235 else
236 {
237 int i, j;
238 int npoints[2] = {m_base[0]->GetNumPoints(), m_base[1]->GetNumPoints()};
239 int nmodes[2] = {m_base[0]->GetNumModes(), m_base[1]->GetNumModes()};
240
241 fill(outarray.get(), outarray.get() + m_ncoeffs, 0.0);
242
243 Array<OneD, NekDouble> physEdge[4];
244 Array<OneD, NekDouble> coeffEdge[4];
245 for (i = 0; i < 4; i++)
246 {
247 physEdge[i] = Array<OneD, NekDouble>(npoints[i % 2]);
248 coeffEdge[i] = Array<OneD, NekDouble>(nmodes[i % 2]);
249 }
250
251 for (i = 0; i < npoints[0]; i++)
252 {
253 physEdge[0][i] = inarray[i];
254 physEdge[2][i] = inarray[npoints[0] * npoints[1] - 1 - i];
255 }
256
257 for (i = 0; i < npoints[1]; i++)
258 {
259 physEdge[1][i] = inarray[npoints[0] - 1 + i * npoints[0]];
260 physEdge[3][i] =
261 inarray[(npoints[1] - 1) * npoints[0] - i * npoints[0]];
262 }
263
264 StdSegExpSharedPtr segexp[2] = {
266 m_base[0]->GetBasisKey()),
268 m_base[1]->GetBasisKey())};
269
270 Array<OneD, unsigned int> mapArray;
271 Array<OneD, int> signArray;
273
274 for (i = 0; i < 4; i++)
275 {
276 segexp[i % 2]->FwdTransBndConstrained(physEdge[i], coeffEdge[i]);
277
278 GetTraceToElementMap(i, mapArray, signArray);
279 for (j = 0; j < nmodes[i % 2]; j++)
280 {
281 sign = (NekDouble)signArray[j];
282 outarray[mapArray[j]] = sign * coeffEdge[i][j];
283 }
284 }
285
286 Array<OneD, NekDouble> tmp0(m_ncoeffs);
287 Array<OneD, NekDouble> tmp1(m_ncoeffs);
288
289 StdMatrixKey masskey(eMass, DetShapeType(), *this);
290 MassMatrixOp(outarray, tmp0, masskey);
291 IProductWRTBase(inarray, tmp1);
292
293 Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
294
295 // get Mass matrix inverse (only of interior DOF)
296 // use block (1,1) of the static condensed system
297 // note: this block alreay contains the inverse matrix
298 DNekMatSharedPtr matsys =
299 (m_stdStaticCondMatrixManager[masskey])->GetBlock(1, 1);
300
301 int nBoundaryDofs = NumBndryCoeffs();
302 int nInteriorDofs = m_ncoeffs - nBoundaryDofs;
303
304 Array<OneD, NekDouble> rhs(nInteriorDofs);
305 Array<OneD, NekDouble> result(nInteriorDofs);
306
307 GetInteriorMap(mapArray);
308
309 for (i = 0; i < nInteriorDofs; i++)
310 {
311 rhs[i] = tmp1[mapArray[i]];
312 }
313
314 Blas::Dgemv('N', nInteriorDofs, nInteriorDofs, 1.0,
315 &(matsys->GetPtr())[0], nInteriorDofs, rhs.get(), 1, 0.0,
316 result.get(), 1);
317
318 for (i = 0; i < nInteriorDofs; i++)
319 {
320 outarray[mapArray[i]] = result[i];
321 }
322 }
323}
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:47
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:752
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:528
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:684
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:674
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 = alpha A x plus beta y where A[m x n].
Definition: Blas.hpp:211
std::shared_ptr< StdSegExp > StdSegExpSharedPtr
Definition: StdSegExp.h:222
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.hpp:220

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 1195 of file StdQuadExp.cpp.

1196{
1197 int i, j;
1198 int order0 = GetBasisNumModes(0);
1199 int order1 = GetBasisNumModes(1);
1200 MatrixType mtype = mkey.GetMatrixType();
1201
1202 DNekMatSharedPtr Mat;
1203
1204 switch (mtype)
1205 {
1207 {
1208 int nq0 = m_base[0]->GetNumPoints();
1209 int nq1 = m_base[1]->GetNumPoints();
1210 int nq;
1211
1212 // take definition from key
1213 if (mkey.ConstFactorExists(eFactorConst))
1214 {
1215 nq = (int)mkey.GetConstFactor(eFactorConst);
1216 }
1217 else
1218 {
1219 nq = max(nq0, nq1);
1220 }
1221
1222 int neq =
1224 Array<OneD, Array<OneD, NekDouble>> coords(neq);
1225 Array<OneD, NekDouble> coll(2);
1226 Array<OneD, DNekMatSharedPtr> I(2);
1227 Array<OneD, NekDouble> tmp(nq0);
1228
1229 Mat = MemoryManager<DNekMat>::AllocateSharedPtr(neq, nq0 * nq1);
1230 int cnt = 0;
1231
1232 for (i = 0; i < nq; ++i)
1233 {
1234 for (j = 0; j < nq; ++j, ++cnt)
1235 {
1236 coords[cnt] = Array<OneD, NekDouble>(2);
1237 coords[cnt][0] = -1.0 + 2 * j / (NekDouble)(nq - 1);
1238 coords[cnt][1] = -1.0 + 2 * i / (NekDouble)(nq - 1);
1239 }
1240 }
1241
1242 for (i = 0; i < neq; ++i)
1243 {
1244 LocCoordToLocCollapsed(coords[i], coll);
1245
1246 I[0] = m_base[0]->GetI(coll);
1247 I[1] = m_base[1]->GetI(coll + 1);
1248
1249 // interpolate first coordinate direction
1250 for (j = 0; j < nq1; ++j)
1251 {
1252 NekDouble fac = (I[1]->GetPtr())[j];
1253 Vmath::Smul(nq0, fac, I[0]->GetPtr(), 1, tmp, 1);
1254
1255 Vmath::Vcopy(nq0, &tmp[0], 1,
1256 Mat->GetRawPtr() + j * nq0 * neq + i, neq);
1257 }
1258 }
1259 break;
1260 }
1261 case eMass:
1262 {
1264 // For Fourier basis set the imaginary component of mean mode
1265 // to have a unit diagonal component in mass matrix
1267 {
1268 for (i = 0; i < order1; ++i)
1269 {
1270 (*Mat)(order0 *i + 1, i * order0 + 1) = 1.0;
1271 }
1272 }
1273
1275 {
1276 for (i = 0; i < order0; ++i)
1277 {
1278 (*Mat)(order0 + i, order0 + i) = 1.0;
1279 }
1280 }
1281 break;
1282 }
1283 case eFwdTrans:
1284 {
1285 Mat =
1287 StdMatrixKey iprodkey(eIProductWRTBase, DetShapeType(), *this);
1288 DNekMat &Iprod = *GetStdMatrix(iprodkey);
1289 StdMatrixKey imasskey(eInvMass, DetShapeType(), *this);
1290 DNekMat &Imass = *GetStdMatrix(imasskey);
1291
1292 (*Mat) = Imass * Iprod;
1293 break;
1294 }
1295 case eGaussDG:
1296 {
1297 ConstFactorMap factors = mkey.GetConstFactors();
1298
1299 int edge = (int)factors[StdRegions::eFactorGaussEdge];
1300 int dir = (edge + 1) % 2;
1301 int nCoeffs = m_base[dir]->GetNumModes();
1302
1303 const LibUtilities::PointsKey BS_p(
1305 const LibUtilities::BasisKey BS_k(LibUtilities::eGauss_Lagrange,
1306 nCoeffs, BS_p);
1307
1308 Array<OneD, NekDouble> coords(1, 0.0);
1309 coords[0] = (edge == 0 || edge == 3) ? -1.0 : 1.0;
1310
1313 DNekMatSharedPtr m_Ix = basis->GetI(coords);
1314
1316 Vmath::Vcopy(nCoeffs, m_Ix->GetPtr(), 1, Mat->GetPtr(), 1);
1317 break;
1318 }
1319 default:
1320 {
1322 break;
1323 }
1324 }
1325
1326 return Mat;
1327}
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:156
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
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:169
int getNumberOfCoefficients(int Na, int Nb)
Definition: ShapeType.hpp:133
BasisManagerT & BasisManager(void)
std::shared_ptr< Basis > BasisSharedPtr
@ eGaussGaussLegendre
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:46
@ eGauss_Lagrange
Lagrange Polynomials using the Gauss points.
Definition: BasisType.h:57
@ eFourier
Fourier Expansion .
Definition: BasisType.h:55
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:430
StdRegions::ConstFactorMap factors
NekMatrix< NekDouble, StandardMatrixTag > DNekMat
Definition: NekTypeDefs.hpp:50
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.hpp:100

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::VarcoeffHashingTest::factors, 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)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 728 of file StdQuadExp.cpp.

729{
730 int i;
731 int cnt = 0;
732 int nummodes0, nummodes1;
733 int value1 = 0, value2 = 0;
734 if (outarray.size() != NumBndryCoeffs())
735 {
736 outarray = Array<OneD, unsigned int>(NumBndryCoeffs());
737 }
738
739 nummodes0 = m_base[0]->GetNumModes();
740 nummodes1 = m_base[1]->GetNumModes();
741
742 const LibUtilities::BasisType Btype0 = GetBasisType(0);
743 const LibUtilities::BasisType Btype1 = GetBasisType(1);
744
745 switch (Btype1)
746 {
749 value1 = nummodes0;
750 break;
752 value1 = 2 * nummodes0;
753 break;
754 default:
755 ASSERTL0(0, "Mapping array is not defined for this expansion");
756 break;
757 }
758
759 for (i = 0; i < value1; i++)
760 {
761 outarray[i] = i;
762 }
763 cnt = value1;
764
765 switch (Btype0)
766 {
769 value2 = value1 + nummodes0 - 1;
770 break;
772 value2 = value1 + 1;
773 break;
774 default:
775 ASSERTL0(0, "Mapping array is not defined for this expansion");
776 break;
777 }
778
779 for (i = 0; i < nummodes1 - 2; i++)
780 {
781 outarray[cnt++] = value1 + i * nummodes0;
782 outarray[cnt++] = value2 + i * nummodes0;
783 }
784
785 if (Btype1 == LibUtilities::eGLL_Lagrange ||
787 {
788 for (i = nummodes0 * (nummodes1 - 1); i < GetNcoeffs(); i++)
789 {
790 outarray[cnt++] = i;
791 }
792 }
793}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
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:56
@ 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 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 680 of file StdQuadExp.cpp.

683{
684 Array<OneD, const NekDouble> z0 = m_base[0]->GetZ();
685 Array<OneD, const NekDouble> z1 = m_base[1]->GetZ();
686 int nq0 = GetNumPoints(0);
687 int nq1 = GetNumPoints(1);
688 int i;
689
690 for (i = 0; i < nq1; ++i)
691 {
692 Blas::Dcopy(nq0, z0.get(), 1, &coords_0[0] + i * nq0, 1);
693 Vmath::Fill(nq0, z1[i], &coords_1[0] + i * nq0, 1);
694 }
695}
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:128
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.hpp:54

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 795 of file StdQuadExp.cpp.

796{
797 int i, j;
798 int cnt = 0;
799 int nummodes0, nummodes1;
800 int startvalue = 0;
801 if (outarray.size() != GetNcoeffs() - NumBndryCoeffs())
802 {
803 outarray = Array<OneD, unsigned int>(GetNcoeffs() - NumBndryCoeffs());
804 }
805
806 nummodes0 = m_base[0]->GetNumModes();
807 nummodes1 = m_base[1]->GetNumModes();
808
809 const LibUtilities::BasisType Btype0 = GetBasisType(0);
810 const LibUtilities::BasisType Btype1 = GetBasisType(1);
811
812 switch (Btype1)
813 {
815 startvalue = nummodes0;
816 break;
818 startvalue = 2 * nummodes0;
819 break;
820 default:
821 ASSERTL0(0, "Mapping array is not defined for this expansion");
822 break;
823 }
824
825 switch (Btype0)
826 {
828 startvalue++;
829 break;
831 startvalue += 2;
832 break;
833 default:
834 ASSERTL0(0, "Mapping array is not defined for this expansion");
835 break;
836 }
837
838 for (i = 0; i < nummodes1 - 2; i++)
839 {
840 for (j = 0; j < nummodes0 - 2; j++)
841 {
842 outarray[cnt++] = startvalue + j;
843 }
844 startvalue += nummodes0;
845 }
846}

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
finalprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 560 of file StdQuadExp.cpp.

561{
562 return 4;
563}

◆ v_GetNverts()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 555 of file StdQuadExp.cpp.

556{
557 return 4;
558}

◆ v_GetSimplexEquiSpacedConnectivity()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1591 of file StdQuadExp.cpp.

1593{
1594 int np1 = m_base[0]->GetNumPoints();
1595 int np2 = m_base[1]->GetNumPoints();
1596 int np = max(np1, np2);
1597
1598 conn = Array<OneD, int>(6 * (np - 1) * (np - 1));
1599
1600 int row = 0;
1601 int rowp1 = 0;
1602 int cnt = 0;
1603 for (int i = 0; i < np - 1; ++i)
1604 {
1605 rowp1 += np;
1606 for (int j = 0; j < np - 1; ++j)
1607 {
1608 conn[cnt++] = row + j;
1609 conn[cnt++] = row + j + 1;
1610 conn[cnt++] = rowp1 + j;
1611
1612 conn[cnt++] = rowp1 + j + 1;
1613 conn[cnt++] = rowp1 + j;
1614 conn[cnt++] = row + j + 1;
1615 }
1616 row += np;
1617 }
1618}

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

◆ v_GetTraceBasisKey()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 606 of file StdQuadExp.cpp.

608{
609 ASSERTL2((i >= 0) && (i <= 3), "edge id is out of range");
610
611 if ((i == 0) || (i == 2))
612 {
613 return GetBasis(0)->GetBasisKey();
614 }
615 else
616 {
617 return GetBasis(1)->GetBasisKey();
618 }
619}
const LibUtilities::BasisSharedPtr & GetBasis(int dir) const
This function gets the shared point to basis in the dir direction.
Definition: StdExpansion.h:112

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

◆ v_GetTraceCoeffMap()

void Nektar::StdRegions::StdQuadExp::v_GetTraceCoeffMap ( const unsigned int  traceid,
Array< OneD, unsigned int > &  maparray 
)
overrideprotectedvirtual

Get the map of the coefficient location to teh local trace coefficients.

Reimplemented from Nektar::StdRegions::StdExpansion2D.

Definition at line 963 of file StdQuadExp.cpp.

965{
966 ASSERTL1(traceid < 4, "traceid must be between 0 and 3");
967
968 unsigned int i;
969 unsigned int order0 = m_base[0]->GetNumModes();
970 unsigned int order1 = m_base[1]->GetNumModes();
971 unsigned int numModes = (traceid % 2) ? order1 : order0;
972
973 if (maparray.size() != numModes)
974 {
975 maparray = Array<OneD, unsigned int>(numModes);
976 }
977
978 const LibUtilities::BasisType bType = GetBasisType(traceid % 2);
979
980 if (bType == LibUtilities::eModified_A)
981 {
982 switch (traceid)
983 {
984 case 0:
985 {
986 for (i = 0; i < numModes; i++)
987 {
988 maparray[i] = i;
989 }
990 }
991 break;
992 case 1:
993 {
994 for (i = 0; i < numModes; i++)
995 {
996 maparray[i] = i * order0 + 1;
997 }
998 }
999 break;
1000 case 2:
1001 {
1002 for (i = 0; i < numModes; i++)
1003 {
1004 maparray[i] = order0 + i;
1005 }
1006 }
1007 break;
1008 case 3:
1009 {
1010 for (i = 0; i < numModes; i++)
1011 {
1012 maparray[i] = i * order0;
1013 }
1014 }
1015 break;
1016 default:
1017 break;
1018 }
1019 }
1020 else if (bType == LibUtilities::eGLL_Lagrange ||
1022 {
1023 switch (traceid)
1024 {
1025 case 0:
1026 {
1027 for (i = 0; i < numModes; i++)
1028 {
1029 maparray[i] = i;
1030 }
1031 }
1032 break;
1033 case 1:
1034 {
1035 for (i = 0; i < numModes; i++)
1036 {
1037 maparray[i] = (i + 1) * order0 - 1;
1038 }
1039 }
1040 break;
1041 case 2:
1042 {
1043 for (i = 0; i < numModes; i++)
1044 {
1045 maparray[i] = order0 * (order1 - 1) + i;
1046 }
1047 }
1048 break;
1049 case 3:
1050 {
1051 for (i = 0; i < numModes; i++)
1052 {
1053 maparray[i] = order0 * i;
1054 }
1055 }
1056 break;
1057 default:
1058 break;
1059 }
1060 }
1061 else
1062 {
1063 ASSERTL0(false, "Mapping not defined for this type of basis");
1064 }
1065}

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

◆ v_GetTraceInteriorToElementMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1067 of file StdQuadExp.cpp.

1070{
1071 int i;
1072 const int nummodes0 = m_base[0]->GetNumModes();
1073 const int nummodes1 = m_base[1]->GetNumModes();
1074 const int nEdgeIntCoeffs = GetTraceNcoeffs(eid) - 2;
1075 const LibUtilities::BasisType bType = GetBasisType(eid % 2);
1076
1077 if (maparray.size() != nEdgeIntCoeffs)
1078 {
1079 maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
1080 }
1081
1082 if (signarray.size() != nEdgeIntCoeffs)
1083 {
1084 signarray = Array<OneD, int>(nEdgeIntCoeffs, 1);
1085 }
1086 else
1087 {
1088 fill(signarray.get(), signarray.get() + nEdgeIntCoeffs, 1);
1089 }
1090
1091 if (bType == LibUtilities::eModified_A)
1092 {
1093 switch (eid)
1094 {
1095 case 0:
1096 {
1097 for (i = 0; i < nEdgeIntCoeffs; i++)
1098 {
1099 maparray[i] = i + 2;
1100 }
1101 }
1102 break;
1103 case 1:
1104 {
1105 for (i = 0; i < nEdgeIntCoeffs; i++)
1106 {
1107 maparray[i] = (i + 2) * nummodes0 + 1;
1108 }
1109 }
1110 break;
1111 case 2:
1112 {
1113 for (i = 0; i < nEdgeIntCoeffs; i++)
1114 {
1115 maparray[i] = nummodes0 + i + 2;
1116 }
1117 }
1118 break;
1119 case 3:
1120 {
1121 for (i = 0; i < nEdgeIntCoeffs; i++)
1122 {
1123 maparray[i] = (i + 2) * nummodes0;
1124 }
1125 }
1126 break;
1127 default:
1128 ASSERTL0(false, "eid must be between 0 and 3");
1129 break;
1130 }
1131
1132 if (edgeOrient == eBackwards)
1133 {
1134 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1135 {
1136 signarray[i] = -1;
1137 }
1138 }
1139 }
1140 else if (bType == LibUtilities::eGLL_Lagrange)
1141 {
1142 switch (eid)
1143 {
1144 case 0:
1145 {
1146 for (i = 0; i < nEdgeIntCoeffs; i++)
1147 {
1148 maparray[i] = i + 1;
1149 }
1150 }
1151 break;
1152 case 1:
1153 {
1154 for (i = 0; i < nEdgeIntCoeffs; i++)
1155 {
1156 maparray[i] = (i + 2) * nummodes0 - 1;
1157 }
1158 }
1159 break;
1160 case 2:
1161 {
1162 for (i = 0; i < nEdgeIntCoeffs; i++)
1163 {
1164 maparray[i] = nummodes0 * (nummodes1 - 1) + i + 1;
1165 }
1166 }
1167 break;
1168 case 3:
1169 {
1170 for (i = 0; i < nEdgeIntCoeffs; i++)
1171 {
1172 maparray[i] = nummodes0 * (i + 1);
1173 }
1174 }
1175 break;
1176 default:
1177 ASSERTL0(false, "eid must be between 0 and 3");
1178 break;
1179 }
1180 if (edgeOrient == eBackwards)
1181 {
1182 reverse(maparray.get(), maparray.get() + nEdgeIntCoeffs);
1183 }
1184 }
1185 else
1186 {
1187 ASSERTL0(false, "Mapping not defined for this type of basis");
1188 }
1189}
int GetTraceNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th trace.
Definition: StdExpansion.h:261

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

int Nektar::StdRegions::StdQuadExp::v_GetTraceIntNcoeffs ( const int  i) const
finalprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 579 of file StdQuadExp.cpp.

580{
581 ASSERTL2((i >= 0) && (i <= 4), "edge id is out of range");
582 if ((i == 0) || (i == 2))
583 {
584 return GetBasisNumModes(0) - 2;
585 }
586 else
587 {
588 return GetBasisNumModes(1) - 2;
589 }
590}

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

◆ v_GetTraceNcoeffs()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 565 of file StdQuadExp.cpp.

566{
567 ASSERTL2((i >= 0) && (i <= 3), "edge id is out of range");
568
569 if ((i == 0) || (i == 2))
570 {
571 return GetBasisNumModes(0);
572 }
573 else
574 {
575 return GetBasisNumModes(1);
576 }
577}

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

◆ v_GetTraceNumPoints()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 592 of file StdQuadExp.cpp.

593{
594 ASSERTL2((i >= 0) && (i <= 3), "edge id is out of range");
595
596 if ((i == 0) || (i == 2))
597 {
598 return GetNumPoints(0);
599 }
600 else
601 {
602 return GetNumPoints(1);
603 }
604}

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

◆ v_GetVertexMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 848 of file StdQuadExp.cpp.

849{
850 int localDOF = 0;
851
852 if (useCoeffPacking == true)
853 {
854 switch (localVertexId)
855 {
856 case 0:
857 {
858 localDOF = 0;
859 }
860 break;
861 case 1:
862 {
864 {
865 localDOF = m_base[0]->GetNumModes() - 1;
866 }
867 else
868 {
869 localDOF = 1;
870 }
871 }
872 break;
873 case 2:
874 {
876 {
877 localDOF = m_base[0]->GetNumModes() *
878 (m_base[1]->GetNumModes() - 1);
879 }
880 else
881 {
882 localDOF = m_base[0]->GetNumModes();
883 }
884 }
885 break;
886 case 3:
887 {
889 {
890 localDOF =
891 m_base[0]->GetNumModes() * m_base[1]->GetNumModes() - 1;
892 }
893 else
894 {
895 localDOF = m_base[0]->GetNumModes() + 1;
896 }
897 }
898 break;
899 default:
900 ASSERTL0(false, "eid must be between 0 and 3");
901 break;
902 }
903 }
904 else
905 {
906 switch (localVertexId)
907 {
908 case 0:
909 {
910 localDOF = 0;
911 }
912 break;
913 case 1:
914 {
916 {
917 localDOF = m_base[0]->GetNumModes() - 1;
918 }
919 else
920 {
921 localDOF = 1;
922 }
923 }
924 break;
925 case 2:
926 {
928 {
929 localDOF =
930 m_base[0]->GetNumModes() * m_base[1]->GetNumModes() - 1;
931 }
932 else
933 {
934 localDOF = m_base[0]->GetNumModes() + 1;
935 }
936 }
937 break;
938 case 3:
939 {
941 {
942 localDOF = m_base[0]->GetNumModes() *
943 (m_base[1]->GetNumModes() - 1);
944 }
945 else
946 {
947 localDOF = m_base[0]->GetNumModes();
948 }
949 }
950 break;
951 default:
952 ASSERTL0(false, "eid must be between 0 and 3");
953 break;
954 }
955 }
956 return localDOF;
957}

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 1558 of file StdQuadExp.cpp.

1561{
1562 StdQuadExp::v_HelmholtzMatrixOp_MatFree(inarray, outarray, mkey);
1563}
void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override

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

◆ v_Integral()

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

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

59{
60 Array<OneD, const NekDouble> w0 = m_base[0]->GetW();
61 Array<OneD, const NekDouble> w1 = m_base[1]->GetW();
62
63 return StdExpansion2D::Integral(inarray, w0, w1);
64}
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 
)
overrideprotectedvirtual

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

356{
357 if (m_base[0]->Collocation() && m_base[1]->Collocation())
358 {
359 MultiplyByQuadratureMetric(inarray, outarray);
360 }
361 else
362 {
363 StdQuadExp::v_IProductWRTBase_SumFac(inarray, outarray);
364 }
365}
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:723
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
Definition: StdQuadExp.cpp:367

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

Referenced by v_FwdTrans().

◆ v_IProductWRTBase_SumFac()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 367 of file StdQuadExp.cpp.

370{
371 int nquad0 = m_base[0]->GetNumPoints();
372 int nquad1 = m_base[1]->GetNumPoints();
373 int order0 = m_base[0]->GetNumModes();
374
375 if (multiplybyweights)
376 {
377 Array<OneD, NekDouble> tmp(nquad0 * nquad1 + nquad1 * order0);
378 Array<OneD, NekDouble> wsp(tmp + nquad0 * nquad1);
379
380 // multiply by integration constants
381 MultiplyByQuadratureMetric(inarray, tmp);
383 m_base[1]->GetBdata(), tmp, outarray, wsp,
384 true, true);
385 }
386 else
387 {
388 Array<OneD, NekDouble> wsp(nquad1 * order0);
390 m_base[1]->GetBdata(), inarray, outarray,
391 wsp, true, true);
392 }
393}
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 
)
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion2D.

Definition at line 443 of file StdQuadExp.cpp.

449{
450 int nquad0 = m_base[0]->GetNumPoints();
451 int nquad1 = m_base[1]->GetNumPoints();
452 int nmodes0 = m_base[0]->GetNumModes();
453 int nmodes1 = m_base[1]->GetNumModes();
454
455 bool colldir0 = doCheckCollDir0 ? (m_base[0]->Collocation()) : false;
456 bool colldir1 = doCheckCollDir1 ? (m_base[1]->Collocation()) : false;
457
458 if (colldir0 && colldir1)
459 {
460 Vmath::Vcopy(m_ncoeffs, inarray.get(), 1, outarray.get(), 1);
461 }
462 else if (colldir0)
463 {
464 Blas::Dgemm('N', 'N', nmodes0, nmodes1, nquad1, 1.0, inarray.get(),
465 nmodes0, base1.get(), nquad1, 0.0, outarray.get(), nmodes0);
466 }
467 else if (colldir1)
468 {
469 Blas::Dgemm('T', 'N', nmodes0, nquad1, nquad0, 1.0, base0.get(), nquad0,
470 inarray.get(), nquad0, 0.0, outarray.get(), nmodes0);
471 }
472 else
473 {
474 ASSERTL1(wsp.size() >= nquad1 * nmodes0,
475 "Workspace size is not sufficient");
476
477#if 1
478 Blas::Dgemm('T', 'N', nmodes0, nquad1, nquad0, 1.0, base0.get(), nquad0,
479 inarray.get(), nquad0, 0.0, wsp.get(), nmodes0);
480
481#else
482 for (int i = 0; i < nmodes0; ++i)
483 {
484 for (int j = 0; j < nquad1; ++j)
485 {
486 wsp[j * nmodes0 + i] =
487 Blas::Ddot(nquad0, base0.get() + i * nquad0, 1,
488 inarray.get() + j * nquad0, 1);
489 }
490 }
491#endif
492 Blas::Dgemm('N', 'N', nmodes0, nmodes1, nquad1, 1.0, wsp.get(), nmodes0,
493 base1.get(), nquad1, 0.0, outarray.get(), nmodes0);
494 }
495}
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:163

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 395 of file StdQuadExp.cpp.

398{
399 v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
400}
void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdQuadExp.cpp:402

References v_IProductWRTDerivBase_SumFac().

◆ v_IProductWRTDerivBase_SumFac()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 402 of file StdQuadExp.cpp.

405{
406 ASSERTL0((dir == 0) || (dir == 1), "input dir is out of range");
407
408 int nquad0 = m_base[0]->GetNumPoints();
409 int nquad1 = m_base[1]->GetNumPoints();
410 int nqtot = nquad0 * nquad1;
411 int order0 = m_base[0]->GetNumModes();
412
413 Array<OneD, NekDouble> tmp(nqtot + nquad1 * order0);
414 Array<OneD, NekDouble> wsp(tmp + nqtot);
415
416 // multiply by integration constants
417 MultiplyByQuadratureMetric(inarray, tmp);
418
419 if (dir) // dir == 1
420 {
422 m_base[1]->GetDbdata(), tmp, outarray, wsp,
423 true, false);
424 }
425 else // dir == 0
426 {
427 IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
428 m_base[1]->GetBdata(), tmp, outarray, wsp,
429 false, true);
430 }
431}

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

Referenced by v_IProductWRTDerivBase().

◆ v_IsBoundaryInteriorExpansion()

bool Nektar::StdRegions::StdQuadExp::v_IsBoundaryInteriorExpansion ( ) const
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 663 of file StdQuadExp.cpp.

664{
665 bool returnval = false;
666
669 {
672 {
673 returnval = true;
674 }
675 }
676
677 return returnval;
678}

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 1537 of file StdQuadExp.cpp.

1540{
1541 StdQuadExp::v_LaplacianMatrixOp_MatFree(inarray, outarray, mkey);
1542}
void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 1544 of file StdQuadExp.cpp.

1547{
1548 StdExpansion::LaplacianMatrixOp_MatFree(k1, k2, inarray, outarray, mkey);
1549}
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 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 508 of file StdQuadExp.cpp.

510{
511 xi[0] = eta[0];
512 xi[1] = eta[1];
513}

◆ v_LocCoordToLocCollapsed()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 501 of file StdQuadExp.cpp.

503{
504 eta[0] = xi[0];
505 eta[1] = xi[1];
506}

◆ v_MassMatrixOp()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 1530 of file StdQuadExp.cpp.

1533{
1534 StdExpansion::MassMatrixOp_MatFree(inarray, outarray, mkey);
1535}
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 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1566 of file StdQuadExp.cpp.

1569{
1570 int i;
1571 int nquad0 = m_base[0]->GetNumPoints();
1572 int nquad1 = m_base[1]->GetNumPoints();
1573
1574 const Array<OneD, const NekDouble> &w0 = m_base[0]->GetW();
1575 const Array<OneD, const NekDouble> &w1 = m_base[1]->GetW();
1576
1577 // multiply by integration constants
1578 for (i = 0; i < nquad1; ++i)
1579 {
1580 Vmath::Vmul(nquad0, inarray.get() + i * nquad0, 1, w0.get(), 1,
1581 outarray.get() + i * nquad0, 1);
1582 }
1583
1584 for (i = 0; i < nquad0; ++i)
1585 {
1586 Vmath::Vmul(nquad1, outarray.get() + i, nquad0, w1.get(), 1,
1587 outarray.get() + i, nquad0);
1588 }
1589}

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

◆ v_NumBndryCoeffs()

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

◆ v_NumDGBndryCoeffs()

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

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

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

80{
81 PhysTensorDeriv(inarray, out_d0, out_d1);
82}
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 
)
overrideprotectedvirtual

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

87{
88 switch (dir)
89 {
90 case 0:
91 {
92 PhysTensorDeriv(inarray, outarray, NullNekDouble1DArray);
93 }
94 break;
95 case 1:
96 {
97 PhysTensorDeriv(inarray, NullNekDouble1DArray, outarray);
98 }
99 break;
100 default:
101 {
102 ASSERTL1(false, "input dir is out of range");
103 }
104 break;
105 }
106}
static Array< OneD, NekDouble > NullNekDouble1DArray

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

◆ v_PhysEvaluate()

NekDouble Nektar::StdRegions::StdQuadExp::v_PhysEvaluate ( const Array< OneD, NekDouble > &  coord,
const Array< OneD, const NekDouble > &  inarray,
std::array< NekDouble, 3 > &  firstOrderDerivs 
)
inlineoverrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion2D.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 155 of file StdQuadExp.h.

159 {
160 return BaryTensorDeriv(coord, inarray, firstOrderDerivs);
161 }
NekDouble BaryTensorDeriv(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)

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

◆ v_PhysEvaluateBasis()

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

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

711{
712 ASSERTL2(coords[0] > -1 - NekConstants::kNekZeroTol, "coord[0] < -1");
713 ASSERTL2(coords[0] < 1 + NekConstants::kNekZeroTol, "coord[0] > 1");
714 ASSERTL2(coords[1] > -1 - NekConstants::kNekZeroTol, "coord[1] < -1");
715 ASSERTL2(coords[1] < 1 + NekConstants::kNekZeroTol, "coord[1] > 1");
716
717 const int nm0 = m_base[0]->GetNumModes();
718 const int nm1 = m_base[1]->GetNumModes();
719
720 return StdExpansion::BaryEvaluateBasis<0>(coords[0], mode % nm1) *
721 StdExpansion::BaryEvaluateBasis<1>(coords[1], mode / nm0);
722}
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 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 1487 of file StdQuadExp.cpp.

1490{
1491 int n_coeffs = inarray.size();
1492
1493 Array<OneD, NekDouble> coeff(n_coeffs);
1494 Array<OneD, NekDouble> coeff_tmp(n_coeffs, 0.0);
1495 Array<OneD, NekDouble> tmp;
1496 Array<OneD, NekDouble> tmp2;
1497
1498 int nmodes0 = m_base[0]->GetNumModes();
1499 int nmodes1 = m_base[1]->GetNumModes();
1500 int numMax = nmodes0;
1501
1502 Vmath::Vcopy(n_coeffs, inarray, 1, coeff_tmp, 1);
1503
1504 const LibUtilities::PointsKey Pkey0(nmodes0,
1506 const LibUtilities::PointsKey Pkey1(nmodes1,
1508
1509 LibUtilities::BasisKey b0(m_base[0]->GetBasisType(), nmodes0, Pkey0);
1510 LibUtilities::BasisKey b1(m_base[1]->GetBasisType(), nmodes1, Pkey1);
1511
1512 LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A, nmodes0, Pkey0);
1513 LibUtilities::BasisKey bortho1(LibUtilities::eOrtho_A, nmodes1, Pkey1);
1514
1515 LibUtilities::InterpCoeff2D(b0, b1, coeff_tmp, bortho0, bortho1, coeff);
1516
1517 Vmath::Zero(n_coeffs, coeff_tmp, 1);
1518
1519 int cnt = 0;
1520 for (int i = 0; i < numMin + 1; ++i)
1521 {
1522 Vmath::Vcopy(numMin, tmp = coeff + cnt, 1, tmp2 = coeff_tmp + cnt, 1);
1523
1524 cnt = i * numMax;
1525 }
1526
1527 LibUtilities::InterpCoeff2D(bortho0, bortho1, coeff_tmp, b0, b1, outarray);
1528}
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:67
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 108 of file StdQuadExp.cpp.

112{
113 StdQuadExp::v_PhysDeriv(inarray, out_d0, out_d1);
114}
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) override
Calculate the derivative of the physical points.
Definition: StdQuadExp.cpp:76

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 116 of file StdQuadExp.cpp.

119{
120 StdQuadExp::v_PhysDeriv(dir, inarray, outarray);
121}

References v_PhysDeriv().

◆ v_SVVLaplacianFilter()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 1338 of file StdQuadExp.cpp.

1340{
1341 // Generate an orthonogal expansion
1342 int qa = m_base[0]->GetNumPoints();
1343 int qb = m_base[1]->GetNumPoints();
1344 int nmodes_a = m_base[0]->GetNumModes();
1345 int nmodes_b = m_base[1]->GetNumModes();
1346 int nmodes = min(nmodes_a, nmodes_b);
1347 // Declare orthogonal basis.
1348 LibUtilities::PointsKey pa(qa, m_base[0]->GetPointsType());
1349 LibUtilities::PointsKey pb(qb, m_base[1]->GetPointsType());
1350
1351 LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A, nmodes_a, pa);
1352 LibUtilities::BasisKey Bb(LibUtilities::eOrtho_A, nmodes_b, pb);
1353 StdQuadExp OrthoExp(Ba, Bb);
1354
1355 // SVV parameters loaded from the .xml case file
1356 Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
1357
1358 // project onto modal space.
1359 OrthoExp.FwdTrans(array, orthocoeffs);
1360
1361 if (mkey.ConstFactorExists(
1362 eFactorSVVPowerKerDiffCoeff)) // Rodrigo's power kernel
1363 {
1364 NekDouble cutoff = mkey.GetConstFactor(eFactorSVVCutoffRatio);
1365 NekDouble SvvDiffCoeff =
1366 mkey.GetConstFactor(eFactorSVVPowerKerDiffCoeff) *
1367 mkey.GetConstFactor(eFactorSVVDiffCoeff);
1368
1369 for (int j = 0; j < nmodes_a; ++j)
1370 {
1371 for (int k = 0; k < nmodes_b; ++k)
1372 {
1373 // linear space but makes high modes very negative
1374 NekDouble fac = std::max(
1375 pow((1.0 * j) / (nmodes_a - 1), cutoff * nmodes_a),
1376 pow((1.0 * k) / (nmodes_b - 1), cutoff * nmodes_b));
1377
1378 orthocoeffs[j * nmodes_b + k] *= SvvDiffCoeff * fac;
1379 }
1380 }
1381 }
1382 else if (mkey.ConstFactorExists(
1383 eFactorSVVDGKerDiffCoeff)) // Rodrigo/mansoor's DG kernel
1384 {
1385 NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDGKerDiffCoeff) *
1386 mkey.GetConstFactor(eFactorSVVDiffCoeff);
1387 int max_ab = max(nmodes_a - kSVVDGFiltermodesmin,
1388 nmodes_b - kSVVDGFiltermodesmin);
1389 max_ab = max(max_ab, 0);
1390 max_ab = min(max_ab, kSVVDGFiltermodesmax - kSVVDGFiltermodesmin);
1391
1392 for (int j = 0; j < nmodes_a; ++j)
1393 {
1394 for (int k = 0; k < nmodes_b; ++k)
1395 {
1396 int maxjk = max(j, k);
1397 maxjk = min(maxjk, kSVVDGFiltermodesmax - 1);
1398
1399 orthocoeffs[j * nmodes_b + k] *=
1400 SvvDiffCoeff * kSVVDGFilter[max_ab][maxjk];
1401 }
1402 }
1403 }
1404 else
1405 {
1406 NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDiffCoeff);
1407 // Exponential Kernel implementation
1408 int cutoff = (int)(mkey.GetConstFactor(eFactorSVVCutoffRatio) *
1409 min(nmodes_a, nmodes_b));
1410
1411 //------"New" Version August 22nd '13--------------------
1412 for (int j = 0; j < nmodes_a; ++j)
1413 {
1414 for (int k = 0; k < nmodes_b; ++k)
1415 {
1416 if (j + k >= cutoff) // to filter out only the "high-modes"
1417 {
1418 orthocoeffs[j * nmodes_b + k] *=
1419 (SvvDiffCoeff *
1420 exp(-(j + k - nmodes) * (j + k - nmodes) /
1421 ((NekDouble)((j + k - cutoff + 1) *
1422 (j + k - cutoff + 1)))));
1423 }
1424 else
1425 {
1426 orthocoeffs[j * nmodes_b + k] *= 0.0;
1427 }
1428 }
1429 }
1430 }
1431
1432 // backward transform to physical space
1433 OrthoExp.BwdTrans(orthocoeffs, array);
1434}
const int kSVVDGFiltermodesmin
Definition: StdRegions.hpp:500
const int kSVVDGFiltermodesmax
Definition: StdRegions.hpp:501
const NekDouble kSVVDGFilter[9][11]
Definition: StdRegions.hpp:503

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 1551 of file StdQuadExp.cpp.

1554{
1555 StdExpansion::WeakDerivMatrixOp_MatFree(i, inarray, outarray, mkey);
1556}
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().