Nektar++
Loading...
Searching...
No Matches
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.
 
 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.
 
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.
 
 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.
 
 StdExpansion (const StdExpansion &T)
 Copy Constructor.
 
virtual ~StdExpansion ()
 Destructor.
 
int GetNumBases () const
 This function returns the number of 1D bases used in the expansion.
 
const Array< OneD, const LibUtilities::BasisSharedPtr > & GetBase () const
 This function gets the shared point to basis.
 
const LibUtilities::BasisSharedPtrGetBasis (int dir) const
 This function gets the shared point to basis in the dir direction.
 
int GetNcoeffs (void) const
 This function returns the total number of coefficients used in the expansion.
 
int GetTotPoints () const
 This function returns the total number of quadrature points used in the element.
 
LibUtilities::BasisType GetBasisType (const int dir) const
 This function returns the type of basis used in the dir direction.
 
int GetBasisNumModes (const int dir) const
 This function returns the number of expansion modes in the dir direction.
 
int EvalBasisNumModesMax (void) const
 This function returns the maximum number of expansion modes over all local directions.
 
LibUtilities::PointsType GetPointsType (const int dir) const
 This function returns the type of quadrature points used in the dir direction.
 
int GetNumPoints (const int dir) const
 This function returns the number of quadrature points in the dir direction.
 
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.
 
int GetNverts () const
 This function returns the number of vertices of the expansion domain.
 
int GetTraceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th trace.
 
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.
 
const LibUtilities::BasisKey GetTraceBasisKey (const int i, int k=-1, bool UseGLL=false) const
 This function returns the basis key belonging to the i-th trace.
 
LibUtilities::PointsKey GetTracePointsKey (const int i, int k=-1) const
 This function returns the basis key belonging to the i-th trace.
 
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.
 
int GetNtraces () const
 Returns the number of trace elements connected to this element.
 
LibUtilities::ShapeType DetShapeType () const
 This function returns the shape of the expansion domain.
 
std::shared_ptr< StdExpansionGetStdExp () const
 
std::shared_ptr< StdExpansionGetLinStdExp (void) const
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion () const
 
bool IsNodalNonTensorialExp ()
 
void NodalToModal (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space.
 
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space.
 
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.
 
void FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 This function fills the array outarray with the mode-th mode of the expansion.
 
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
 
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.
 
void SetElmtId (const int id)
 Set the element id of this expansion when used in a list by returning value of m_elmt_id.
 
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
 
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
 
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}\)
 
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)
 
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.
 
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.
 
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.
 
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.
 
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta.
 
void LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi)
 Convert local collapsed coordinates eta into local cartesian coordinate xi.
 
void PhysInterp (std::shared_ptr< StdExpansion > fromExp, const Array< OneD, const NekDouble > &fromData, Array< OneD, NekDouble > &toData)
 interpolate from one set of quadrature points available from FromExp to the set of quadrature points in the current expansion. If the points are the same this routine will just copy the data
 
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.
 
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.
 
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.
 
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.
 
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.
 
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.
 
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.
 
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.
 
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.
 
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_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.
 
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.
 
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.
 
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, bool UseGLL=false) 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.
 
NekDouble v_PhysEvalFirstDeriv (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.
 
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.
 
NekDouble v_PhysEvaluateInterp (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals) override
 
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_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.
 
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
 
void v_PhysInterp (std::shared_ptr< StdExpansion > fromExp, const Array< OneD, const NekDouble > &fromData, Array< OneD, NekDouble > &toData) 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.
 
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 NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
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.
 
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.
 
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 120 of file StdQuadExp.cpp.

122{
123 if (m_base[0]->Collocation() && m_base[1]->Collocation())
124 {
126 inarray, 1, outarray, 1);
127 }
128 else
129 {
130 StdQuadExp::v_BwdTrans_SumFac(inarray, outarray);
131 }
132}
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
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 134 of file StdQuadExp.cpp.

136{
137 Array<OneD, NekDouble> wsp(m_base[0]->GetNumPoints() *
138 m_base[1]->GetNumModes());
139
140 BwdTrans_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetBdata(), inarray,
141 outarray, wsp, true, true);
142}
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 154 of file StdQuadExp.cpp.

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

679{
680 int nmodes = nummodes[modes_offset] * nummodes[modes_offset + 1];
681 modes_offset += 2;
682
683 return nmodes;
684}

◆ v_CreateStdMatrix()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1360 of file StdQuadExp.cpp.

1361{
1362 return GenMatrix(mkey);
1363}
DNekMatSharedPtr GenMatrix(const StdMatrixKey &mkey)

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

1471{
1472 // Generate an orthogonal expansion
1473 int qa = m_base[0]->GetNumPoints();
1474 int qb = m_base[1]->GetNumPoints();
1475 int nmodesA = m_base[0]->GetNumModes();
1476 int nmodesB = m_base[1]->GetNumModes();
1477 int P = nmodesA - 1;
1478 int Q = nmodesB - 1;
1479
1480 // Declare orthogonal basis.
1481 LibUtilities::PointsKey pa(qa, m_base[0]->GetPointsType());
1482 LibUtilities::PointsKey pb(qb, m_base[1]->GetPointsType());
1483
1484 LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A, nmodesA, pa);
1485 LibUtilities::BasisKey Bb(LibUtilities::eOrtho_A, nmodesB, pb);
1486 StdQuadExp OrthoExp(Ba, Bb);
1487
1488 // Cutoff
1489 int Pcut = cutoff * P;
1490 int Qcut = cutoff * Q;
1491
1492 // Project onto orthogonal space.
1493 Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
1494 OrthoExp.FwdTrans(array, orthocoeffs);
1495
1496 //
1497 NekDouble fac, fac1, fac2;
1498 for (int i = 0; i < nmodesA; ++i)
1499 {
1500 for (int j = 0; j < nmodesB; ++j)
1501 {
1502 // to filter out only the "high-modes"
1503 if (i > Pcut || j > Qcut)
1504 {
1505 fac1 = (NekDouble)(i - Pcut) / ((NekDouble)(P - Pcut));
1506 fac2 = (NekDouble)(j - Qcut) / ((NekDouble)(Q - Qcut));
1507 fac = max(fac1, fac2);
1508 fac = pow(fac, exponent);
1509 orthocoeffs[i * nmodesB + j] *= exp(-alpha * fac);
1510 }
1511 }
1512 }
1513
1514 // backward transform to physical space
1515 OrthoExp.BwdTrans(orthocoeffs, array);
1516}
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
StdQuadExp(const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb)
Constructor using BasisKey class for quadrature points and order definition.
@ P
Monomial polynomials .
Definition BasisType.h:62
@ eOrtho_A
Principle Orthogonal Functions .
Definition BasisType.h:42
scalarT< T > max(scalarT< T > lhs, scalarT< T > rhs)
Definition scalar.hpp:305

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, tinysimd::max(), 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 517 of file StdQuadExp.cpp.

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

Definition at line 197 of file StdQuadExp.cpp.

199{
200 if ((m_base[0]->Collocation()) && (m_base[1]->Collocation()))
201 {
202 Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
203 }
204 else
205 {
206 StdQuadExp::v_IProductWRTBase(inarray, outarray);
207
208 // get Mass matrix inverse
209 StdMatrixKey masskey(eInvMass, DetShapeType(), *this);
210 DNekMatSharedPtr matsys = GetStdMatrix(masskey);
211
212 // copy inarray in case inarray == outarray
213 NekVector<NekDouble> in(m_ncoeffs, outarray, eCopy);
214 NekVector<NekDouble> out(m_ncoeffs, outarray, eWrapper);
215
216 out = (*matsys) * in;
217 }
218}
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
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.
std::shared_ptr< DNekMat > DNekMatSharedPtr

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.

Definition at line 220 of file StdQuadExp.cpp.

223{
224 if ((m_base[0]->Collocation()) && (m_base[1]->Collocation()))
225 {
226 Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
227 }
228 else
229 {
230 int i, j;
231 int npoints[2] = {m_base[0]->GetNumPoints(), m_base[1]->GetNumPoints()};
232 int nmodes[2] = {m_base[0]->GetNumModes(), m_base[1]->GetNumModes()};
233
234 fill(outarray.data(), outarray.data() + m_ncoeffs, 0.0);
235
236 Array<OneD, NekDouble> physEdge[4];
237 Array<OneD, NekDouble> coeffEdge[4];
238 for (i = 0; i < 4; i++)
239 {
240 physEdge[i] = Array<OneD, NekDouble>(npoints[i % 2]);
241 coeffEdge[i] = Array<OneD, NekDouble>(nmodes[i % 2]);
242 }
243
244 for (i = 0; i < npoints[0]; i++)
245 {
246 physEdge[0][i] = inarray[i];
247 physEdge[2][i] = inarray[npoints[0] * npoints[1] - 1 - i];
248 }
249
250 for (i = 0; i < npoints[1]; i++)
251 {
252 physEdge[1][i] = inarray[npoints[0] - 1 + i * npoints[0]];
253 physEdge[3][i] =
254 inarray[(npoints[1] - 1) * npoints[0] - i * npoints[0]];
255 }
256
257 StdSegExpSharedPtr segexp[2] = {
259 m_base[0]->GetBasisKey()),
261 m_base[1]->GetBasisKey())};
262
263 Array<OneD, unsigned int> mapArray;
264 Array<OneD, int> signArray;
266
267 for (i = 0; i < 4; i++)
268 {
269 segexp[i % 2]->FwdTransBndConstrained(physEdge[i], coeffEdge[i]);
270
271 GetTraceToElementMap(i, mapArray, signArray);
272 for (j = 0; j < nmodes[i % 2]; j++)
273 {
274 sign = (NekDouble)signArray[j];
275 outarray[mapArray[j]] = sign * coeffEdge[i][j];
276 }
277 }
278
279 Array<OneD, NekDouble> tmp0(m_ncoeffs);
280 Array<OneD, NekDouble> tmp1(m_ncoeffs);
281
282 StdMatrixKey masskey(eMass, DetShapeType(), *this);
283 MassMatrixOp(outarray, tmp0, masskey);
284 IProductWRTBase(inarray, tmp1);
285
286 Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
287
288 // get Mass matrix inverse (only of interior DOF)
289 // use block (1,1) of the static condensed system
290 // note: this block alreay contains the inverse matrix
291 DNekMatSharedPtr matsys =
292 (m_stdStaticCondMatrixManager[masskey])->GetBlock(1, 1);
293
294 int nBoundaryDofs = NumBndryCoeffs();
295 int nInteriorDofs = m_ncoeffs - nBoundaryDofs;
296
297 Array<OneD, NekDouble> rhs(nInteriorDofs);
298 Array<OneD, NekDouble> result(nInteriorDofs);
299
300 GetInteriorMap(mapArray);
301
302 for (i = 0; i < nInteriorDofs; i++)
303 {
304 rhs[i] = tmp1[mapArray[i]];
305 }
306
307 Blas::Dgemv('N', nInteriorDofs, nInteriorDofs, 1.0,
308 &(matsys->GetPtr())[0], nInteriorDofs, rhs.data(), 1, 0.0,
309 result.data(), 1);
310
311 for (i = 0; i < nInteriorDofs; i++)
312 {
313 outarray[mapArray[i]] = result[i];
314 }
315 }
316}
#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)
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...
void GetTraceToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
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:214
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.

Definition at line 1226 of file StdQuadExp.cpp.

1227{
1228 int i, j;
1229 int order0 = GetBasisNumModes(0);
1230 int order1 = GetBasisNumModes(1);
1231 MatrixType mtype = mkey.GetMatrixType();
1232
1233 DNekMatSharedPtr Mat;
1234
1235 switch (mtype)
1236 {
1238 {
1239 int nq0 = m_base[0]->GetNumPoints();
1240 int nq1 = m_base[1]->GetNumPoints();
1241 int nq;
1242
1243 // take definition from key
1244 if (mkey.ConstFactorExists(eFactorConst))
1245 {
1246 nq = (int)mkey.GetConstFactor(eFactorConst);
1247 }
1248 else
1249 {
1250 nq = max(nq0, nq1);
1251 }
1252
1253 int neq =
1255 Array<OneD, Array<OneD, NekDouble>> coords(neq);
1256 Array<OneD, NekDouble> coll(2);
1257 Array<OneD, DNekMatSharedPtr> I(2);
1258 Array<OneD, NekDouble> tmp(nq0);
1259
1260 Mat = MemoryManager<DNekMat>::AllocateSharedPtr(neq, nq0 * nq1);
1261 int cnt = 0;
1262
1263 for (i = 0; i < nq; ++i)
1264 {
1265 for (j = 0; j < nq; ++j, ++cnt)
1266 {
1267 coords[cnt] = Array<OneD, NekDouble>(2);
1268 coords[cnt][0] = -1.0 + 2 * j / (NekDouble)(nq - 1);
1269 coords[cnt][1] = -1.0 + 2 * i / (NekDouble)(nq - 1);
1270 }
1271 }
1272
1273 for (i = 0; i < neq; ++i)
1274 {
1275 LocCoordToLocCollapsed(coords[i], coll);
1276
1277 I[0] = m_base[0]->GetI(coll);
1278 I[1] = m_base[1]->GetI(coll + 1);
1279
1280 // interpolate first coordinate direction
1281 for (j = 0; j < nq1; ++j)
1282 {
1283 NekDouble fac = (I[1]->GetPtr())[j];
1284 Vmath::Smul(nq0, fac, I[0]->GetPtr(), 1, tmp, 1);
1285
1286 Vmath::Vcopy(nq0, &tmp[0], 1,
1287 Mat->GetRawPtr() + j * nq0 * neq + i, neq);
1288 }
1289 }
1290 break;
1291 }
1292 case eMass:
1293 {
1295 // For Fourier basis set the imaginary component of mean mode
1296 // to have a unit diagonal component in mass matrix
1298 {
1299 for (i = 0; i < order1; ++i)
1300 {
1301 (*Mat)(order0 *i + 1, i * order0 + 1) = 1.0;
1302 }
1303 }
1304
1306 {
1307 for (i = 0; i < order0; ++i)
1308 {
1309 (*Mat)(order0 + i, order0 + i) = 1.0;
1310 }
1311 }
1312 break;
1313 }
1314 case eFwdTrans:
1315 {
1316 Mat =
1318 StdMatrixKey iprodkey(eIProductWRTBase, DetShapeType(), *this);
1319 DNekMat &Iprod = *GetStdMatrix(iprodkey);
1320 StdMatrixKey imasskey(eInvMass, DetShapeType(), *this);
1321 DNekMat &Imass = *GetStdMatrix(imasskey);
1322
1323 (*Mat) = Imass * Iprod;
1324 break;
1325 }
1326 case eGaussDG:
1327 {
1328 ConstFactorMap factors = mkey.GetConstFactors();
1329
1330 int edge = (int)factors[StdRegions::eFactorGaussEdge];
1331 int dir = (edge + 1) % 2;
1332 int nCoeffs = m_base[dir]->GetNumModes();
1333
1334 const LibUtilities::PointsKey BS_p(
1336 const LibUtilities::BasisKey BS_k(LibUtilities::eGauss_Lagrange,
1337 nCoeffs, BS_p);
1338
1339 Array<OneD, NekDouble> coords(1, 0.0);
1340 coords[0] = (edge == 0 || edge == 3) ? -1.0 : 1.0;
1341
1344 DNekMatSharedPtr m_Ix = basis->GetI(coords);
1345
1347 Vmath::Vcopy(nCoeffs, m_Ix->GetPtr(), 1, Mat->GetPtr(), 1);
1348 break;
1349 }
1350 default:
1351 {
1353 break;
1354 }
1355 }
1356
1357 return Mat;
1358}
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
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.
int getNumberOfCoefficients(int Na, int Nb)
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
StdRegions::ConstFactorMap factors
NekMatrix< NekDouble, StandardMatrixTag > DNekMat
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::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, tinysimd::max(), 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 759 of file StdQuadExp.cpp.

760{
761 int i;
762 int cnt = 0;
763 int nummodes0, nummodes1;
764 int value1 = 0, value2 = 0;
765 if (outarray.size() != NumBndryCoeffs())
766 {
767 outarray = Array<OneD, unsigned int>(NumBndryCoeffs());
768 }
769
770 nummodes0 = m_base[0]->GetNumModes();
771 nummodes1 = m_base[1]->GetNumModes();
772
773 const LibUtilities::BasisType Btype0 = GetBasisType(0);
774 const LibUtilities::BasisType Btype1 = GetBasisType(1);
775
776 switch (Btype1)
777 {
780 value1 = nummodes0;
781 break;
783 value1 = 2 * nummodes0;
784 break;
785 default:
786 ASSERTL0(0, "Mapping array is not defined for this expansion");
787 break;
788 }
789
790 for (i = 0; i < value1; i++)
791 {
792 outarray[i] = i;
793 }
794 cnt = value1;
795
796 switch (Btype0)
797 {
800 value2 = value1 + nummodes0 - 1;
801 break;
803 value2 = value1 + 1;
804 break;
805 default:
806 ASSERTL0(0, "Mapping array is not defined for this expansion");
807 break;
808 }
809
810 for (i = 0; i < nummodes1 - 2; i++)
811 {
812 outarray[cnt++] = value1 + i * nummodes0;
813 outarray[cnt++] = value2 + i * nummodes0;
814 }
815
816 if (Btype1 == LibUtilities::eGLL_Lagrange ||
818 {
819 for (i = nummodes0 * (nummodes1 - 1); i < GetNcoeffs(); i++)
820 {
821 outarray[cnt++] = i;
822 }
823 }
824}
#define ASSERTL0(condition, msg)
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
@ 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.

Definition at line 703 of file StdQuadExp.cpp.

706{
707 Array<OneD, const NekDouble> z0 = m_base[0]->GetZ();
708 Array<OneD, const NekDouble> z1 = m_base[1]->GetZ();
709 int nq0 = GetNumPoints(0);
710 int nq1 = GetNumPoints(1);
711 int i;
712
713 for (i = 0; i < nq1; ++i)
714 {
715 Blas::Dcopy(nq0, z0.data(), 1, &coords_0[0] + i * nq0, 1);
716 Vmath::Fill(nq0, z1[i], &coords_1[0] + i * nq0, 1);
717 }
718}
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 826 of file StdQuadExp.cpp.

827{
828 int i, j;
829 int cnt = 0;
830 int nummodes0, nummodes1;
831 int startvalue = 0;
832 if (outarray.size() != GetNcoeffs() - NumBndryCoeffs())
833 {
834 outarray = Array<OneD, unsigned int>(GetNcoeffs() - NumBndryCoeffs());
835 }
836
837 nummodes0 = m_base[0]->GetNumModes();
838 nummodes1 = m_base[1]->GetNumModes();
839
840 const LibUtilities::BasisType Btype0 = GetBasisType(0);
841 const LibUtilities::BasisType Btype1 = GetBasisType(1);
842
843 switch (Btype1)
844 {
846 startvalue = nummodes0;
847 break;
849 startvalue = 2 * nummodes0;
850 break;
851 default:
852 ASSERTL0(0, "Mapping array is not defined for this expansion");
853 break;
854 }
855
856 switch (Btype0)
857 {
859 startvalue++;
860 break;
862 startvalue += 2;
863 break;
864 default:
865 ASSERTL0(0, "Mapping array is not defined for this expansion");
866 break;
867 }
868
869 for (i = 0; i < nummodes1 - 2; i++)
870 {
871 for (j = 0; j < nummodes0 - 2; j++)
872 {
873 outarray[cnt++] = startvalue + j;
874 }
875 startvalue += nummodes0;
876 }
877}

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

557{
558 return 4;
559}

◆ v_GetNverts()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 551 of file StdQuadExp.cpp.

552{
553 return 4;
554}

◆ v_GetSimplexEquiSpacedConnectivity()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1622 of file StdQuadExp.cpp.

1624{
1625 int np1 = m_base[0]->GetNumPoints();
1626 int np2 = m_base[1]->GetNumPoints();
1627 int np = max(np1, np2);
1628
1629 conn = Array<OneD, int>(6 * (np - 1) * (np - 1));
1630
1631 int row = 0;
1632 int rowp1 = 0;
1633 int cnt = 0;
1634 for (int i = 0; i < np - 1; ++i)
1635 {
1636 rowp1 += np;
1637 for (int j = 0; j < np - 1; ++j)
1638 {
1639 conn[cnt++] = row + j;
1640 conn[cnt++] = row + j + 1;
1641 conn[cnt++] = rowp1 + j;
1642
1643 conn[cnt++] = rowp1 + j + 1;
1644 conn[cnt++] = rowp1 + j;
1645 conn[cnt++] = row + j + 1;
1646 }
1647 row += np;
1648 }
1649}

References Nektar::StdRegions::StdExpansion::m_base, and tinysimd::max().

◆ v_GetTraceBasisKey()

const LibUtilities::BasisKey Nektar::StdRegions::StdQuadExp::v_GetTraceBasisKey ( const int  i,
const int  j,
bool  UseGLL = false 
) const
finalprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 602 of file StdQuadExp.cpp.

605{
606 ASSERTL2((i >= 0) && (i <= 3), "edge id is out of range");
607
608 if ((i == 0) || (i == 2))
609 {
610 switch (GetBasis(0)->GetBasisType())
611 {
613 {
614 return LibUtilities::BasisKey(LibUtilities::eGLL_Lagrange,
615 GetBasis(0)->GetNumModes(),
616 GetBasis(0)->GetPointsKey());
617 }
618 break;
619 default:
620 {
621 return GetBasis(0)->GetBasisKey();
622 }
623 }
624 }
625 else
626 {
627 switch (GetBasis(1)->GetBasisType())
628 {
630 {
631 return LibUtilities::BasisKey(LibUtilities::eGLL_Lagrange,
632 GetBasis(1)->GetNumModes(),
633 GetBasis(1)->GetPointsKey());
634 }
635 break;
636 default:
637 {
638 return GetBasis(1)->GetBasisKey();
639 }
640 }
641 }
642}
const LibUtilities::BasisSharedPtr & GetBasis(int dir) const
This function gets the shared point to basis in the dir direction.

References ASSERTL2, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eOrtho_A, Nektar::StdRegions::StdExpansion::GetBasis(), and Nektar::StdRegions::StdExpansion::GetBasisType().

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

996{
997 ASSERTL1(traceid < 4, "traceid must be between 0 and 3");
998
999 unsigned int i;
1000 unsigned int order0 = m_base[0]->GetNumModes();
1001 unsigned int order1 = m_base[1]->GetNumModes();
1002 unsigned int numModes = (traceid % 2) ? order1 : order0;
1003
1004 if (maparray.size() != numModes)
1005 {
1006 maparray = Array<OneD, unsigned int>(numModes);
1007 }
1008
1009 const LibUtilities::BasisType bType = GetBasisType(traceid % 2);
1010
1011 if (bType == LibUtilities::eModified_A)
1012 {
1013 switch (traceid)
1014 {
1015 case 0:
1016 {
1017 for (i = 0; i < numModes; i++)
1018 {
1019 maparray[i] = i;
1020 }
1021 }
1022 break;
1023 case 1:
1024 {
1025 for (i = 0; i < numModes; i++)
1026 {
1027 maparray[i] = i * order0 + 1;
1028 }
1029 }
1030 break;
1031 case 2:
1032 {
1033 for (i = 0; i < numModes; i++)
1034 {
1035 maparray[i] = order0 + i;
1036 }
1037 }
1038 break;
1039 case 3:
1040 {
1041 for (i = 0; i < numModes; i++)
1042 {
1043 maparray[i] = i * order0;
1044 }
1045 }
1046 break;
1047 default:
1048 break;
1049 }
1050 }
1051 else if (bType == LibUtilities::eGLL_Lagrange ||
1053 {
1054 switch (traceid)
1055 {
1056 case 0:
1057 {
1058 for (i = 0; i < numModes; i++)
1059 {
1060 maparray[i] = i;
1061 }
1062 }
1063 break;
1064 case 1:
1065 {
1066 for (i = 0; i < numModes; i++)
1067 {
1068 maparray[i] = (i + 1) * order0 - 1;
1069 }
1070 }
1071 break;
1072 case 2:
1073 {
1074 for (i = 0; i < numModes; i++)
1075 {
1076 maparray[i] = order0 * (order1 - 1) + i;
1077 }
1078 }
1079 break;
1080 case 3:
1081 {
1082 for (i = 0; i < numModes; i++)
1083 {
1084 maparray[i] = order0 * i;
1085 }
1086 }
1087 break;
1088 default:
1089 break;
1090 }
1091 }
1092 else
1093 {
1094 ASSERTL0(false, "Mapping not defined for this type of basis");
1095 }
1096}

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

1101{
1102 int i;
1103 const int nummodes0 = m_base[0]->GetNumModes();
1104 const int nummodes1 = m_base[1]->GetNumModes();
1105 const int nEdgeIntCoeffs = GetTraceNcoeffs(eid) - 2;
1106 const LibUtilities::BasisType bType = GetBasisType(eid % 2);
1107
1108 if (maparray.size() != nEdgeIntCoeffs)
1109 {
1110 maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
1111 }
1112
1113 if (signarray.size() != nEdgeIntCoeffs)
1114 {
1115 signarray = Array<OneD, int>(nEdgeIntCoeffs, 1);
1116 }
1117 else
1118 {
1119 fill(signarray.data(), signarray.data() + nEdgeIntCoeffs, 1);
1120 }
1121
1122 if (bType == LibUtilities::eModified_A)
1123 {
1124 switch (eid)
1125 {
1126 case 0:
1127 {
1128 for (i = 0; i < nEdgeIntCoeffs; i++)
1129 {
1130 maparray[i] = i + 2;
1131 }
1132 }
1133 break;
1134 case 1:
1135 {
1136 for (i = 0; i < nEdgeIntCoeffs; i++)
1137 {
1138 maparray[i] = (i + 2) * nummodes0 + 1;
1139 }
1140 }
1141 break;
1142 case 2:
1143 {
1144 for (i = 0; i < nEdgeIntCoeffs; i++)
1145 {
1146 maparray[i] = nummodes0 + i + 2;
1147 }
1148 }
1149 break;
1150 case 3:
1151 {
1152 for (i = 0; i < nEdgeIntCoeffs; i++)
1153 {
1154 maparray[i] = (i + 2) * nummodes0;
1155 }
1156 }
1157 break;
1158 default:
1159 ASSERTL0(false, "eid must be between 0 and 3");
1160 break;
1161 }
1162
1163 if (edgeOrient == eBackwards)
1164 {
1165 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1166 {
1167 signarray[i] = -1;
1168 }
1169 }
1170 }
1171 else if (bType == LibUtilities::eGLL_Lagrange)
1172 {
1173 switch (eid)
1174 {
1175 case 0:
1176 {
1177 for (i = 0; i < nEdgeIntCoeffs; i++)
1178 {
1179 maparray[i] = i + 1;
1180 }
1181 }
1182 break;
1183 case 1:
1184 {
1185 for (i = 0; i < nEdgeIntCoeffs; i++)
1186 {
1187 maparray[i] = (i + 2) * nummodes0 - 1;
1188 }
1189 }
1190 break;
1191 case 2:
1192 {
1193 for (i = 0; i < nEdgeIntCoeffs; i++)
1194 {
1195 maparray[i] = nummodes0 * (nummodes1 - 1) + i + 1;
1196 }
1197 }
1198 break;
1199 case 3:
1200 {
1201 for (i = 0; i < nEdgeIntCoeffs; i++)
1202 {
1203 maparray[i] = nummodes0 * (i + 1);
1204 }
1205 }
1206 break;
1207 default:
1208 ASSERTL0(false, "eid must be between 0 and 3");
1209 break;
1210 }
1211 if (edgeOrient == eBackwards)
1212 {
1213 reverse(maparray.data(), maparray.data() + nEdgeIntCoeffs);
1214 }
1215 }
1216 else
1217 {
1218 ASSERTL0(false, "Mapping not defined for this type of basis");
1219 }
1220}
int GetTraceNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th trace.

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

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

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

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

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

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

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

880{
881 int localDOF = 0;
882
883 if (useCoeffPacking == true)
884 {
885 switch (localVertexId)
886 {
887 case 0:
888 {
889 localDOF = 0;
890 }
891 break;
892 case 1:
893 {
895 {
896 localDOF = m_base[0]->GetNumModes() - 1;
897 }
898 else
899 {
900 localDOF = 1;
901 }
902 }
903 break;
904 case 2:
905 {
907 {
908 localDOF = m_base[0]->GetNumModes() *
909 (m_base[1]->GetNumModes() - 1);
910 }
911 else
912 {
913 localDOF = m_base[0]->GetNumModes();
914 }
915 }
916 break;
917 case 3:
918 {
920 {
921 localDOF =
922 m_base[0]->GetNumModes() * m_base[1]->GetNumModes() - 1;
923 }
924 else
925 {
926 localDOF = m_base[0]->GetNumModes() + 1;
927 }
928 }
929 break;
930 default:
931 ASSERTL0(false, "eid must be between 0 and 3");
932 break;
933 }
934 }
935 else
936 {
937 switch (localVertexId)
938 {
939 case 0:
940 {
941 localDOF = 0;
942 }
943 break;
944 case 1:
945 {
947 {
948 localDOF = m_base[0]->GetNumModes() - 1;
949 }
950 else
951 {
952 localDOF = 1;
953 }
954 }
955 break;
956 case 2:
957 {
959 {
960 localDOF =
961 m_base[0]->GetNumModes() * m_base[1]->GetNumModes() - 1;
962 }
963 else
964 {
965 localDOF = m_base[0]->GetNumModes() + 1;
966 }
967 }
968 break;
969 case 3:
970 {
972 {
973 localDOF = m_base[0]->GetNumModes() *
974 (m_base[1]->GetNumModes() - 1);
975 }
976 else
977 {
978 localDOF = m_base[0]->GetNumModes();
979 }
980 }
981 break;
982 default:
983 ASSERTL0(false, "eid must be between 0 and 3");
984 break;
985 }
986 }
987 return localDOF;
988}

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.

Definition at line 1589 of file StdQuadExp.cpp.

1592{
1593 StdQuadExp::v_HelmholtzMatrixOp_MatFree(inarray, outarray, mkey);
1594}
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)

References Nektar::StdRegions::StdExpansion::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.

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.

Definition at line 347 of file StdQuadExp.cpp.

349{
350 if (m_base[0]->Collocation() && m_base[1]->Collocation())
351 {
352 MultiplyByQuadratureMetric(inarray, outarray);
353 }
354 else
355 {
356 StdQuadExp::v_IProductWRTBase_SumFac(inarray, outarray);
357 }
358}
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override

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.

Definition at line 360 of file StdQuadExp.cpp.

363{
364 int nquad0 = m_base[0]->GetNumPoints();
365 int nquad1 = m_base[1]->GetNumPoints();
366 int order0 = m_base[0]->GetNumModes();
367
368 if (multiplybyweights)
369 {
370 Array<OneD, NekDouble> tmp(nquad0 * nquad1 + nquad1 * order0);
371 Array<OneD, NekDouble> wsp(tmp + nquad0 * nquad1);
372
373 // multiply by integration constants
374 MultiplyByQuadratureMetric(inarray, tmp);
376 m_base[1]->GetBdata(), tmp, outarray, wsp,
377 true, true);
378 }
379 else
380 {
381 Array<OneD, NekDouble> wsp(nquad1 * order0);
383 m_base[1]->GetBdata(), inarray, outarray,
384 wsp, true, true);
385 }
386}
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 436 of file StdQuadExp.cpp.

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

Definition at line 388 of file StdQuadExp.cpp.

391{
392 v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
393}
void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override

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.

Definition at line 395 of file StdQuadExp.cpp.

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

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

687{
688 bool returnval = false;
689
692 {
695 {
696 returnval = true;
697 }
698 }
699
700 return returnval;
701}

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.

Definition at line 1568 of file StdQuadExp.cpp.

1571{
1572 StdQuadExp::v_LaplacianMatrixOp_MatFree(inarray, outarray, mkey);
1573}
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)

References Nektar::StdRegions::StdExpansion::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.

Definition at line 1575 of file StdQuadExp.cpp.

1578{
1579 StdExpansion::LaplacianMatrixOp_MatFree(k1, k2, inarray, outarray, mkey);
1580}
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 504 of file StdQuadExp.cpp.

506{
507 xi[0] = eta[0];
508 xi[1] = eta[1];
509}

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

499{
500 eta[0] = xi[0];
501 eta[1] = xi[1];
502}

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

Definition at line 1561 of file StdQuadExp.cpp.

1564{
1565 StdExpansion::MassMatrixOp_MatFree(inarray, outarray, mkey);
1566}
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 1597 of file StdQuadExp.cpp.

1600{
1601 int i;
1602 int nquad0 = m_base[0]->GetNumPoints();
1603 int nquad1 = m_base[1]->GetNumPoints();
1604
1605 const Array<OneD, const NekDouble> &w0 = m_base[0]->GetW();
1606 const Array<OneD, const NekDouble> &w1 = m_base[1]->GetW();
1607
1608 // multiply by integration constants
1609 for (i = 0; i < nquad1; ++i)
1610 {
1611 Vmath::Vmul(nquad0, inarray.data() + i * nquad0, 1, w0.data(), 1,
1612 outarray.data() + i * nquad0, 1);
1613 }
1614
1615 for (i = 0; i < nquad0; ++i)
1616 {
1617 Vmath::Vmul(nquad1, outarray.data() + i, nquad0, w1.data(), 1,
1618 outarray.data() + i, nquad0);
1619 }
1620}

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.

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.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 747 of file StdQuadExp.cpp.

751{
752 return BaryTensorDeriv(coord, inarray, firstOrderDerivs);
753}
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 732 of file StdQuadExp.cpp.

734{
735 ASSERTL2(coords[0] > -1 - NekConstants::kNekZeroTol, "coord[0] < -1");
736 ASSERTL2(coords[0] < 1 + NekConstants::kNekZeroTol, "coord[0] > 1");
737 ASSERTL2(coords[1] > -1 - NekConstants::kNekZeroTol, "coord[1] < -1");
738 ASSERTL2(coords[1] < 1 + NekConstants::kNekZeroTol, "coord[1] > 1");
739
740 const int nm0 = m_base[0]->GetNumModes();
741 const int nm1 = m_base[1]->GetNumModes();
742
743 return StdExpansion::BaryEvaluateBasis<0>(coords[0], mode % nm1) *
744 StdExpansion::BaryEvaluateBasis<1>(coords[1], mode / nm0);
745}
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.

Definition at line 1518 of file StdQuadExp.cpp.

1521{
1522 int n_coeffs = inarray.size();
1523
1524 Array<OneD, NekDouble> coeff(n_coeffs);
1525 Array<OneD, NekDouble> coeff_tmp(n_coeffs, 0.0);
1526 Array<OneD, NekDouble> tmp;
1527 Array<OneD, NekDouble> tmp2;
1528
1529 int nmodes0 = m_base[0]->GetNumModes();
1530 int nmodes1 = m_base[1]->GetNumModes();
1531 int numMax = nmodes0;
1532
1533 Vmath::Vcopy(n_coeffs, inarray, 1, coeff_tmp, 1);
1534
1535 const LibUtilities::PointsKey Pkey0(nmodes0,
1537 const LibUtilities::PointsKey Pkey1(nmodes1,
1539
1540 LibUtilities::BasisKey b0(m_base[0]->GetBasisType(), nmodes0, Pkey0);
1541 LibUtilities::BasisKey b1(m_base[1]->GetBasisType(), nmodes1, Pkey1);
1542
1543 LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A, nmodes0, Pkey0);
1544 LibUtilities::BasisKey bortho1(LibUtilities::eOrtho_A, nmodes1, Pkey1);
1545
1546 LibUtilities::InterpCoeff2D(b0, b1, coeff_tmp, bortho0, bortho1, coeff);
1547
1548 Vmath::Zero(n_coeffs, coeff_tmp, 1);
1549
1550 int cnt = 0;
1551 for (int i = 0; i < numMin + 1; ++i)
1552 {
1553 Vmath::Vcopy(numMin, tmp = coeff + cnt, 1, tmp2 = coeff_tmp + cnt, 1);
1554
1555 cnt = i * numMax;
1556 }
1557
1558 LibUtilities::InterpCoeff2D(bortho0, bortho1, coeff_tmp, b0, b1, outarray);
1559}
void InterpCoeff2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
@ 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()

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.

References v_PhysDeriv().

◆ v_SVVLaplacianFilter()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1369 of file StdQuadExp.cpp.

1371{
1372 // Generate an orthonogal expansion
1373 int qa = m_base[0]->GetNumPoints();
1374 int qb = m_base[1]->GetNumPoints();
1375 int nmodes_a = m_base[0]->GetNumModes();
1376 int nmodes_b = m_base[1]->GetNumModes();
1377 int nmodes = min(nmodes_a, nmodes_b);
1378 // Declare orthogonal basis.
1379 LibUtilities::PointsKey pa(qa, m_base[0]->GetPointsType());
1380 LibUtilities::PointsKey pb(qb, m_base[1]->GetPointsType());
1381
1382 LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A, nmodes_a, pa);
1383 LibUtilities::BasisKey Bb(LibUtilities::eOrtho_A, nmodes_b, pb);
1384 StdQuadExp OrthoExp(Ba, Bb);
1385
1386 // SVV parameters loaded from the .xml case file
1387 Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
1388
1389 // project onto modal space.
1390 OrthoExp.FwdTrans(array, orthocoeffs);
1391
1392 if (mkey.ConstFactorExists(
1393 eFactorSVVPowerKerDiffCoeff)) // Rodrigo's power kernel
1394 {
1395 NekDouble cutoff = mkey.GetConstFactor(eFactorSVVCutoffRatio);
1396 NekDouble SvvDiffCoeff =
1397 mkey.GetConstFactor(eFactorSVVPowerKerDiffCoeff) *
1398 mkey.GetConstFactor(eFactorSVVDiffCoeff);
1399
1400 for (int j = 0; j < nmodes_a; ++j)
1401 {
1402 for (int k = 0; k < nmodes_b; ++k)
1403 {
1404 // linear space but makes high modes very negative
1405 NekDouble fac = std::max(
1406 pow((1.0 * j) / (nmodes_a - 1), cutoff * nmodes_a),
1407 pow((1.0 * k) / (nmodes_b - 1), cutoff * nmodes_b));
1408
1409 orthocoeffs[j * nmodes_b + k] *= SvvDiffCoeff * fac;
1410 }
1411 }
1412 }
1413 else if (mkey.ConstFactorExists(
1414 eFactorSVVDGKerDiffCoeff)) // Rodrigo/mansoor's DG kernel
1415 {
1416 NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDGKerDiffCoeff) *
1417 mkey.GetConstFactor(eFactorSVVDiffCoeff);
1418 int max_ab = max(nmodes_a - kSVVDGFiltermodesmin,
1419 nmodes_b - kSVVDGFiltermodesmin);
1420 max_ab = max(max_ab, 0);
1421 max_ab = min(max_ab, kSVVDGFiltermodesmax - kSVVDGFiltermodesmin);
1422
1423 for (int j = 0; j < nmodes_a; ++j)
1424 {
1425 for (int k = 0; k < nmodes_b; ++k)
1426 {
1427 int maxjk = max(j, k);
1428 maxjk = min(maxjk, kSVVDGFiltermodesmax - 1);
1429
1430 orthocoeffs[j * nmodes_b + k] *=
1431 SvvDiffCoeff * kSVVDGFilter[max_ab][maxjk];
1432 }
1433 }
1434 }
1435 else
1436 {
1437 NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDiffCoeff);
1438 // Exponential Kernel implementation
1439 int cutoff = (int)(mkey.GetConstFactor(eFactorSVVCutoffRatio) *
1440 min(nmodes_a, nmodes_b));
1441
1442 //------"New" Version August 22nd '13--------------------
1443 for (int j = 0; j < nmodes_a; ++j)
1444 {
1445 for (int k = 0; k < nmodes_b; ++k)
1446 {
1447 if (j + k >= cutoff) // to filter out only the "high-modes"
1448 {
1449 orthocoeffs[j * nmodes_b + k] *=
1450 (SvvDiffCoeff *
1451 exp(-(j + k - nmodes) * (j + k - nmodes) /
1452 ((NekDouble)((j + k - cutoff + 1) *
1453 (j + k - cutoff + 1)))));
1454 }
1455 else
1456 {
1457 orthocoeffs[j * nmodes_b + k] *= 0.0;
1458 }
1459 }
1460 }
1461 }
1462
1463 // backward transform to physical space
1464 OrthoExp.BwdTrans(orthocoeffs, array);
1465}
const int kSVVDGFiltermodesmin
const int kSVVDGFiltermodesmax
const NekDouble kSVVDGFilter[9][11]
scalarT< T > min(scalarT< T > lhs, scalarT< T > rhs)
Definition scalar.hpp:300

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, Nektar::StdRegions::StdExpansion::m_base, tinysimd::max(), and tinysimd::min().

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

Definition at line 1582 of file StdQuadExp.cpp.

1585{
1586 StdExpansion::WeakDerivMatrixOp_MatFree(i, inarray, outarray, mkey);
1587}
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().