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

Class representing a segment element in reference space All interface of this class sits in StdExpansion class. More...

#include <StdSegExp.h>

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

Public Member Functions

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

Protected Member Functions

virtual NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray) override
 Integrate the physical point list inarray over region and return the value. More...
 
void v_PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray) override
 Evaluate the derivative \( d/d{\xi_1} \) at the physical quadrature points given by inarray and return in outarray. More...
 
void v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Calculate the derivative of the physical points in a given direction. More...
 
void v_StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray) override
 
void v_StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Backward transform from coefficient space given in inarray and evaluate at the physical quadrature points outarray. More...
 
void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Forward transform from physical quadrature space stored in inarray and evaluate the expansion coefficients and store in outarray. More...
 
void v_BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
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
 Inner product of inarray over region with respect to the expansion basis (this)->m_base[0] and return in outarray. More...
 
void v_IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check) override
 Inner product of inarray over region with respect to expansion basis base and return in outarray. More...
 
void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) 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
 
virtual void v_LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta) override
 
virtual void v_LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi) override
 
NekDouble v_PhysEvaluateBasis (const Array< OneD, const NekDouble > &coords, int mode) final override
 
NekDouble v_PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
 
NekDouble v_PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs, std::array< NekDouble, 6 > &secondOrderDerivs) override
 
virtual void v_LaplacianMatrixOp (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
 
virtual void v_SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey) override
 
virtual void v_ExponentialFilter (Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff) override
 
virtual void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
virtual void v_FillMode (const int mode, Array< OneD, NekDouble > &outarray) override
 
virtual void v_GetCoords (Array< OneD, NekDouble > &coords_0, Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2) override
 
virtual void v_GetBoundaryMap (Array< OneD, unsigned int > &outarray) override
 
virtual void v_GetInteriorMap (Array< OneD, unsigned int > &outarray) override
 
virtual int v_GetVertexMap (int localVertexId, bool useCoeffPacking=false) override
 
virtual int v_GetNverts () const final override
 
virtual int v_GetNtraces () const final override
 
virtual int v_GetTraceNcoeffs (const int i) const final override
 
virtual int v_GetTraceIntNcoeffs (const int i) const final override
 
virtual int v_GetTraceNumPoints (const int i) const final override
 
virtual int v_NumBndryCoeffs () const override
 
virtual int v_NumDGBndryCoeffs () const override
 
virtual bool v_IsBoundaryInteriorExpansion () const override
 
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset) override
 
virtual LibUtilities::ShapeType v_DetShapeType () const override
 Return Shape of region, using ShapeType enum list. i.e. Segment. More...
 
virtual DNekMatSharedPtr v_GenMatrix (const StdMatrixKey &mkey) override
 
virtual DNekMatSharedPtr v_CreateStdMatrix (const StdMatrixKey &mkey) override
 
virtual void v_GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true) override
 
void v_ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
virtual void v_GetTraceCoeffMap (const unsigned int traceid, Array< OneD, unsigned int > &maparray) override
 Get the map of the coefficient location to teh local trace coefficients. More...
 
void v_GetElmtTraceToTraceMap (const unsigned int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient, int P, int Q) override
 
void v_GetTraceToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient, int P, int Q) override
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion1D
virtual NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals) override
 
virtual NekDouble v_PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
 
virtual NekDouble v_PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs, std::array< NekDouble, 6 > &secondOrderDerivs) override
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion
DNekMatSharedPtr CreateStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr CreateStdStaticCondMatrix (const StdMatrixKey &mkey)
 Create the static condensation of a matrix when using a boundary interior decomposition. More...
 
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase_SumFac (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GeneralMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
 
void LaplacianMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp_MatFree (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void HelmholtzMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
 
virtual NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
virtual void v_GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 
virtual void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals, NekDouble &deriv, NekDouble &deriv2)
 This function performs the barycentric interpolation of the polynomial stored in coord at a point physvals using barycentric interpolation weights in direction. More...
 
template<int DIR>
NekDouble BaryEvaluateBasis (const NekDouble &coord, const int &mode)
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals)
 Helper function to pass an unused value by reference into BaryEvaluate. More...
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals, NekDouble &deriv)
 

Additional Inherited Members

- Protected Attributes inherited from Nektar::StdRegions::StdExpansion
Array< OneD, LibUtilities::BasisSharedPtrm_base
 
int m_elmt_id
 
int m_ncoeffs
 
LibUtilities::NekManager< StdMatrixKey, DNekMat, StdMatrixKey::opLessm_stdMatrixManager
 
LibUtilities::NekManager< StdMatrixKey, DNekBlkMat, StdMatrixKey::opLessm_stdStaticCondMatrixManager
 

Detailed Description

Class representing a segment element in reference space All interface of this class sits in StdExpansion class.

Definition at line 47 of file StdSegExp.h.

Constructor & Destructor Documentation

◆ StdSegExp() [1/3]

Nektar::StdRegions::StdSegExp::StdSegExp ( )

defult constructor

Definition at line 47 of file StdSegExp.cpp.

48{
49}

◆ StdSegExp() [2/3]

Nektar::StdRegions::StdSegExp::StdSegExp ( const LibUtilities::BasisKey Ba)

Constructor using BasisKey class for quadrature points and order definition.

Parameters
BaBasisKey class definition containing order and quadrature points.

Definition at line 57 of file StdSegExp.cpp.

58 : StdExpansion(Ba.GetNumModes(), 1, Ba),
59 StdExpansion1D(Ba.GetNumModes(), Ba)
60{
61}
StdExpansion()
Default Constructor.

◆ StdSegExp() [3/3]

Nektar::StdRegions::StdSegExp::StdSegExp ( const StdSegExp T)

Copy Constructor.

Definition at line 64 of file StdSegExp.cpp.

65{
66}

◆ ~StdSegExp()

Nektar::StdRegions::StdSegExp::~StdSegExp ( )
overridevirtual

Destructor.

Definition at line 69 of file StdSegExp.cpp.

70{
71}

Member Function Documentation

◆ v_BwdTrans()

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

Backward transform from coefficient space given in inarray and evaluate at the physical quadrature points outarray.

Operation can be evaluated as \( u(\xi_{1i}) = \sum_{p=0}^{order-1} \hat{u}_p \phi_p(\xi_{1i}) \) or equivalently \( {\bf u} = {\bf B}^T {\bf \hat{u}} \) where \({\bf B}[i][j] = \phi_i(\xi_{1j}), \mbox{\_coeffs}[p] = {\bf \hat{u}}[p] \)

The function takes the coefficient array inarray as input for the transformation

Parameters
inarraythe coeffficients of the expansion
outarraythe resulting array of the values of the function at the physical quadrature points will be stored in the array outarray

Implements Nektar::StdRegions::StdExpansion.

Definition at line 202 of file StdSegExp.cpp.

204{
205 int nquad = m_base[0]->GetNumPoints();
206
207 if (m_base[0]->Collocation())
208 {
209 Vmath::Vcopy(nquad, inarray, 1, outarray, 1);
210 }
211 else
212 {
213
214 Blas::Dgemv('N', nquad, m_base[0]->GetNumModes(), 1.0,
215 (m_base[0]->GetBdata()).get(), nquad, &inarray[0], 1, 0.0,
216 &outarray[0], 1);
217 }
218}
Array< OneD, LibUtilities::BasisSharedPtr > m_base
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:213
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191

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

Referenced by v_BwdTrans_SumFac(), v_HelmholtzMatrixOp(), and v_LaplacianMatrixOp().

◆ v_BwdTrans_SumFac()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 365 of file StdSegExp.cpp.

367{
368 v_BwdTrans(inarray, outarray);
369}
void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Backward transform from coefficient space given in inarray and evaluate at the physical quadrature po...
Definition: StdSegExp.cpp:202

References v_BwdTrans().

◆ v_CalcNumberOfCoefficients()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 694 of file StdSegExp.cpp.

696{
697 int nmodes = nummodes[modes_offset];
698 modes_offset += 1;
699
700 return nmodes;
701}

◆ v_CreateStdMatrix()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 770 of file StdSegExp.cpp.

771{
772 return v_GenMatrix(mkey);
773}
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override
Definition: StdSegExp.cpp:707

References v_GenMatrix().

◆ v_DetShapeType()

LibUtilities::ShapeType Nektar::StdRegions::StdSegExp::v_DetShapeType ( ) const
overrideprotectedvirtual

Return Shape of region, using ShapeType enum list. i.e. Segment.

Implements Nektar::StdRegions::StdExpansion.

Definition at line 76 of file StdSegExp.cpp.

References Nektar::LibUtilities::eSegment.

Referenced by v_FwdTrans(), v_FwdTransBndConstrained(), and v_GenMatrix().

◆ v_ExponentialFilter()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 593 of file StdSegExp.cpp.

597{
598 // Generate an orthogonal expansion
599 int nq = m_base[0]->GetNumPoints();
600 int nmodes = m_base[0]->GetNumModes();
601 int P = nmodes - 1;
602 // Declare orthogonal basis.
603 LibUtilities::PointsKey pKey(nq, m_base[0]->GetPointsType());
604
605 LibUtilities::BasisKey B(LibUtilities::eOrtho_A, nmodes, pKey);
606 StdSegExp OrthoExp(B);
607
608 // Cutoff
609 int Pcut = cutoff * P;
610
611 // Project onto orthogonal space.
612 Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
613 OrthoExp.FwdTrans(array, orthocoeffs);
614
615 //
616 NekDouble fac;
617 for (int j = 0; j < nmodes; ++j)
618 {
619 // to filter out only the "high-modes"
620 if (j > Pcut)
621 {
622 fac = (NekDouble)(j - Pcut) / ((NekDouble)(P - Pcut));
623 fac = pow(fac, exponent);
624 orthocoeffs[j] *= exp(-alpha * fac);
625 }
626 }
627
628 // backward transform to physical space
629 OrthoExp.BwdTrans(orthocoeffs, array);
630}
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:211
StdSegExp()
defult constructor
Definition: StdSegExp.cpp:47
@ P
Monomial polynomials .
Definition: BasisType.h:64
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:44
double NekDouble

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

◆ v_FillMode()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 495 of file StdSegExp.cpp.

496{
497 int nquad = m_base[0]->GetNumPoints();
498 const NekDouble *base = m_base[0]->GetBdata().get();
499
500 ASSERTL2(mode <= m_ncoeffs,
501 "calling argument mode is larger than total expansion order");
502
503 Vmath::Vcopy(nquad, (NekDouble *)base + mode * nquad, 1, &outarray[0], 1);
504}
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:272

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

◆ v_FwdTrans()

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

Forward transform from physical quadrature space stored in inarray and evaluate the expansion coefficients and store in outarray.

Perform a forward transform using a Galerkin projection by taking the inner product of the physical points and multiplying by the inverse of the mass matrix using the Solve method of the standard matrix container holding the local mass matrix, i.e. \( {\bf \hat{u}} = {\bf M}^{-1} {\bf I} \) where \( {\bf I}[p] = \int^1_{-1} \phi_p(\xi_1) u(\xi_1) d\xi_1 \)

This function stores the expansion coefficients calculated by the transformation in the coefficient space array outarray

Parameters
inarrayarray of physical quadrature points to be transformed
outarraythe coeffficients of the expansion

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 269 of file StdSegExp.cpp.

271{
272 if (m_base[0]->Collocation())
273 {
274 Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
275 }
276 else
277 {
278 v_IProductWRTBase(inarray, outarray);
279
280 // get Mass matrix inverse
281 StdMatrixKey masskey(eInvMass, v_DetShapeType(), *this);
282 DNekMatSharedPtr matsys = GetStdMatrix(masskey);
283
284 NekVector<NekDouble> in(m_ncoeffs, outarray, eCopy);
285 NekVector<NekDouble> out(m_ncoeffs, outarray, eWrapper);
286
287 out = (*matsys) * in;
288 }
289}
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:609
void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Inner product of inarray over region with respect to the expansion basis (this)->m_base[0] and return...
Definition: StdSegExp.cpp:435
virtual LibUtilities::ShapeType v_DetShapeType() const override
Return Shape of region, using ShapeType enum list. i.e. Segment.
Definition: StdSegExp.cpp:76
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75

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

Referenced by v_FwdTransBndConstrained().

◆ v_FwdTransBndConstrained()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 291 of file StdSegExp.cpp.

294{
295 if (m_base[0]->Collocation())
296 {
297 Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
298 }
299 else
300 {
301 int nInteriorDofs = m_ncoeffs - 2;
302 int offset = 0;
303
304 switch (m_base[0]->GetBasisType())
305 {
307 {
308 offset = 1;
309 }
310 break;
312 {
313 nInteriorDofs = m_ncoeffs;
314 offset = 0;
315 }
316 break;
319 {
320 offset = 2;
321 }
322 break;
323 default:
324 ASSERTL0(false, "This type of FwdTrans is not defined for this "
325 "expansion type");
326 }
327
328 fill(outarray.get(), outarray.get() + m_ncoeffs, 0.0);
329
331 {
332 outarray[GetVertexMap(0)] = inarray[0];
333 outarray[GetVertexMap(1)] = inarray[m_base[0]->GetNumPoints() - 1];
334
335 if (m_ncoeffs > 2)
336 {
337 // ideally, we would like to have tmp0 to be replaced by
338 // outarray (currently MassMatrixOp does not allow aliasing)
339 Array<OneD, NekDouble> tmp0(m_ncoeffs);
340 Array<OneD, NekDouble> tmp1(m_ncoeffs);
341
342 StdMatrixKey masskey(eMass, v_DetShapeType(), *this);
343 MassMatrixOp(outarray, tmp0, masskey);
344 v_IProductWRTBase(inarray, tmp1);
345
346 Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
347
348 // get Mass matrix inverse (only of interior DOF)
349 DNekMatSharedPtr matsys =
350 (m_stdStaticCondMatrixManager[masskey])->GetBlock(1, 1);
351
352 Blas::Dgemv('N', nInteriorDofs, nInteriorDofs, 1.0,
353 &(matsys->GetPtr())[0], nInteriorDofs,
354 tmp1.get() + offset, 1, 0.0,
355 outarray.get() + offset, 1);
356 }
357 }
358 else
359 {
360 StdSegExp::v_FwdTrans(inarray, outarray);
361 }
362 }
363}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:162
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:758
int GetVertexMap(const int localVertexId, bool useCoeffPacking=false)
Definition: StdExpansion.h:685
LibUtilities::NekManager< StdMatrixKey, DNekBlkMat, StdMatrixKey::opLess > m_stdStaticCondMatrixManager
void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Forward transform from physical quadrature space stored in inarray and evaluate the expansion coeffic...
Definition: StdSegExp.cpp:269
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:51
@ eGauss_Lagrange
Lagrange Polynomials using the Gauss points.
Definition: BasisType.h:59
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:58
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:50
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:414

References ASSERTL0, Blas::Dgemv(), Nektar::LibUtilities::eGauss_Lagrange, Nektar::LibUtilities::eGLL_Lagrange, Nektar::StdRegions::eMass, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::GetVertexMap(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::StdRegions::StdExpansion::m_stdStaticCondMatrixManager, Nektar::StdRegions::StdExpansion::MassMatrixOp(), v_DetShapeType(), v_FwdTrans(), v_IProductWRTBase(), Vmath::Vcopy(), and Vmath::Vsub().

◆ v_GenMatrix()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 707 of file StdSegExp.cpp.

708{
710 MatrixType mattype;
711
712 switch (mattype = mkey.GetMatrixType())
713 {
715 {
716 int nq = m_base[0]->GetNumPoints();
717
718 // take definition from key
719 if (mkey.ConstFactorExists(eFactorConst))
720 {
721 nq = (int)mkey.GetConstFactor(eFactorConst);
722 }
723
725 Array<OneD, NekDouble> coords(1);
728
729 for (int i = 0; i < neq; ++i)
730 {
731 coords[0] = -1.0 + 2 * i / (NekDouble)(neq - 1);
732 I = m_base[0]->GetI(coords);
733 Vmath::Vcopy(nq, I->GetRawPtr(), 1, Mat->GetRawPtr() + i, neq);
734 }
735 }
736 break;
737 case eFwdTrans:
738 {
739 Mat =
741 StdMatrixKey iprodkey(eIProductWRTBase, v_DetShapeType(), *this);
742 DNekMat &Iprod = *GetStdMatrix(iprodkey);
743 StdMatrixKey imasskey(eInvMass, v_DetShapeType(), *this);
744 DNekMat &Imass = *GetStdMatrix(imasskey);
745
746 (*Mat) = Imass * Iprod;
747 }
748 break;
749 default:
750 {
752
753 if (mattype == eMass)
754 {
755 // For Fourier basis set the imaginary component
756 // of mean mode to have a unit diagonal component
757 // in mass matrix
759 {
760 (*Mat)(1, 1) = 1.0;
761 }
762 }
763 }
764 break;
765 }
766
767 return Mat;
768}
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
@ eFourier
Fourier Expansion .
Definition: BasisType.h:57
NekMatrix< NekDouble, StandardMatrixTag > DNekMat
Definition: NekTypeDefs.hpp:50

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::StdRegions::StdMatrixKey::ConstFactorExists(), Nektar::StdRegions::StdExpansion::CreateGeneralMatrix(), Nektar::StdRegions::eFactorConst, Nektar::LibUtilities::eFourier, Nektar::StdRegions::eFwdTrans, Nektar::StdRegions::eInvMass, Nektar::StdRegions::eIProductWRTBase, Nektar::StdRegions::eMass, Nektar::StdRegions::ePhysInterpToEquiSpaced, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdMatrixKey::GetMatrixType(), Nektar::LibUtilities::StdSegData::getNumberOfCoefficients(), Nektar::StdRegions::StdExpansion::GetStdMatrix(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, v_DetShapeType(), and Vmath::Vcopy().

Referenced by v_CreateStdMatrix().

◆ v_GetBoundaryMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 779 of file StdSegExp.cpp.

780{
781 if (outarray.size() != NumBndryCoeffs())
782 {
783 outarray = Array<OneD, unsigned int>(NumBndryCoeffs());
784 }
785 const LibUtilities::BasisType Btype = GetBasisType(0);
786 int nummodes = m_base[0]->GetNumModes();
787
788 outarray[0] = 0;
789
790 switch (Btype)
791 {
796 outarray[1] = nummodes - 1;
797 break;
800 outarray[1] = 1;
801 break;
802 default:
803 ASSERTL0(0, "Mapping array is not defined for this expansion");
804 break;
805 }
806}
@ eChebyshev
Chebyshev Polynomials.
Definition: BasisType.h:63

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

◆ v_GetCoords()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 644 of file StdSegExp.cpp.

647{
648 boost::ignore_unused(coords_1, coords_2);
649 Blas::Dcopy(GetNumPoints(0), (m_base[0]->GetZ()).get(), 1, &coords_0[0], 1);
650}
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:224
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:130

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

◆ v_GetElmtTraceToTraceMap()

void Nektar::StdRegions::StdSegExp::v_GetElmtTraceToTraceMap ( const unsigned int  eid,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray,
Orientation  edgeOrient,
int  P,
int  Q 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 925 of file StdSegExp.cpp.

929{
930 // parameters for higher dimnesion traces
931 boost::ignore_unused(eid, orient, P, Q);
932 if (maparray.size() != 1)
933 {
934 maparray = Array<OneD, unsigned int>(1);
935 }
936
937 maparray[0] = 0;
938
939 if (signarray.size() != 1)
940 {
941 signarray = Array<OneD, int>(1, 1);
942 }
943 else
944 {
945 signarray[0] = 1;
946 }
947}

References Nektar::LibUtilities::P.

◆ v_GetInteriorMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 808 of file StdSegExp.cpp.

809{
810 int i;
811 if (outarray.size() != GetNcoeffs() - NumBndryCoeffs())
812 {
813 outarray = Array<OneD, unsigned int>(GetNcoeffs() - NumBndryCoeffs());
814 }
815 const LibUtilities::BasisType Btype = GetBasisType(0);
816
817 switch (Btype)
818 {
823 for (i = 0; i < GetNcoeffs() - 2; i++)
824 {
825 outarray[i] = i + 1;
826 }
827 break;
830 for (i = 0; i < GetNcoeffs() - 2; i++)
831 {
832 outarray[i] = i + 2;
833 }
834 break;
835 default:
836 ASSERTL0(0, "Mapping array is not defined for this expansion");
837 break;
838 }
839}
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:130

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

◆ v_GetNtraces()

int Nektar::StdRegions::StdSegExp::v_GetNtraces ( ) const
finaloverrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 661 of file StdSegExp.cpp.

662{
663 return 2;
664}

◆ v_GetNverts()

int Nektar::StdRegions::StdSegExp::v_GetNverts ( ) const
finaloverrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 656 of file StdSegExp.cpp.

657{
658 return 2;
659}

◆ v_GetSimplexEquiSpacedConnectivity()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 858 of file StdSegExp.cpp.

860{
861 boost::ignore_unused(standard);
862 int np = m_base[0]->GetNumPoints();
863
864 conn = Array<OneD, int>(2 * (np - 1));
865 int cnt = 0;
866 for (int i = 0; i < np - 1; ++i)
867 {
868 conn[cnt++] = i;
869 conn[cnt++] = i + 1;
870 }
871}

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

◆ v_GetTraceCoeffMap()

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

Definition at line 877 of file StdSegExp.cpp.

879{
880 int order0 = m_base[0]->GetNumModes();
881
882 ASSERTL0(traceid < 2, "eid must be between 0 and 1");
883
884 if (maparray.size() != 1)
885 {
886 maparray = Array<OneD, unsigned int>(1);
887 }
888
889 const LibUtilities::BasisType bType = GetBasisType(0);
890
891 if (bType == LibUtilities::eModified_A)
892 {
893 maparray[0] = (traceid == 0) ? 0 : 1;
894 }
895 else if (bType == LibUtilities::eGLL_Lagrange ||
897 {
898 maparray[0] = (traceid == 0) ? 0 : order0 - 1;
899 }
900 else
901 {
902 ASSERTL0(false, "Unknown Basis");
903 }
904}

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

Referenced by v_GetTraceToElementMap().

◆ v_GetTraceIntNcoeffs()

int Nektar::StdRegions::StdSegExp::v_GetTraceIntNcoeffs ( const int  i) const
finaloverrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 672 of file StdSegExp.cpp.

673{
674 boost::ignore_unused(i);
675 return 0;
676}

◆ v_GetTraceNcoeffs()

int Nektar::StdRegions::StdSegExp::v_GetTraceNcoeffs ( const int  i) const
finaloverrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 666 of file StdSegExp.cpp.

667{
668 boost::ignore_unused(i);
669 return 1;
670}

◆ v_GetTraceNumPoints()

int Nektar::StdRegions::StdSegExp::v_GetTraceNumPoints ( const int  i) const
finaloverrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 678 of file StdSegExp.cpp.

679{
680 boost::ignore_unused(i);
681 return 1;
682}

◆ v_GetTraceToElementMap()

void Nektar::StdRegions::StdSegExp::v_GetTraceToElementMap ( const int  tid,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray,
Orientation  edgeOrient,
int  P,
int  Q 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 906 of file StdSegExp.cpp.

910{
911 boost::ignore_unused(orient, P, Q);
912
913 v_GetTraceCoeffMap(tid, maparray);
914
915 if (signarray.size() != 1)
916 {
917 signarray = Array<OneD, int>(1, 1);
918 }
919 else
920 {
921 signarray[0] = 1;
922 }
923}
virtual 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.
Definition: StdSegExp.cpp:877

References Nektar::LibUtilities::P, and v_GetTraceCoeffMap().

◆ v_GetVertexMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 841 of file StdSegExp.cpp.

842{
843 boost::ignore_unused(useCoeffPacking);
844 ASSERTL0((localVertexId == 0) || (localVertexId == 1),
845 "local vertex id"
846 "must be between 0 or 1");
847
848 int localDOF = localVertexId;
849
851 (localVertexId == 1))
852 {
853 localDOF = m_base[0]->GetNumModes() - 1;
854 }
855 return localDOF;
856}

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

◆ v_HelmholtzMatrixOp()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 530 of file StdSegExp.cpp.

533{
534 int nquad = m_base[0]->GetNumPoints();
535
536 Array<OneD, NekDouble> physValues(nquad);
537 Array<OneD, NekDouble> dPhysValuesdx(nquad);
538 Array<OneD, NekDouble> wsp(m_ncoeffs);
539
540 v_BwdTrans(inarray, physValues);
541
542 // mass matrix operation
543 v_IProductWRTBase((m_base[0]->GetBdata()), physValues, wsp, 1);
544
545 // Laplacian matrix operation
546 v_PhysDeriv(physValues, dPhysValuesdx);
547 v_IProductWRTBase(m_base[0]->GetDbdata(), dPhysValuesdx, outarray, 1);
548 Blas::Daxpy(m_ncoeffs, mkey.GetConstFactor(eFactorLambda), wsp.get(), 1,
549 outarray.get(), 1);
550}
void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray) override
Evaluate the derivative at the physical quadrature points given by inarray and return in outarray.
Definition: StdSegExp.cpp:140
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
Definition: Blas.hpp:137

References Blas::Daxpy(), Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, v_BwdTrans(), v_IProductWRTBase(), and v_PhysDeriv().

◆ v_Integral()

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

Integrate the physical point list inarray over region and return the value.

Parameters
inarraydefinition of function to be integrated evauluated at quadrature point of expansion.
Returns
returns \(\int^1_{-1} u(\xi_1)d \xi_1 \) where \(inarray[i] = u(\xi_{1i}) \)

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 111 of file StdSegExp.cpp.

112{
113 NekDouble Int = 0.0;
114 int nquad0 = m_base[0]->GetNumPoints();
115 Array<OneD, NekDouble> tmp(nquad0);
116 Array<OneD, const NekDouble> z = m_base[0]->GetZ();
117 Array<OneD, const NekDouble> w0 = m_base[0]->GetW();
118
119 // multiply by integration constants
120 Vmath::Vmul(nquad0, inarray, 1, w0, 1, tmp, 1);
121
122 Int = Vmath::Vsum(nquad0, tmp, 1);
123
124 return Int;
125}
std::vector< double > z(NPUPPER)
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:207
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:890

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

◆ v_IProductWRTBase() [1/2]

void Nektar::StdRegions::StdSegExp::v_IProductWRTBase ( const Array< OneD, const NekDouble > &  base,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
int  coll_check 
)
overrideprotectedvirtual

Inner product of inarray over region with respect to expansion basis base and return in outarray.

Calculate \( I[p] = \int^{1}_{-1} \phi_p(\xi_1) u(\xi_1) d\xi_1 = \sum_{i=0}^{nq-1} \phi_p(\xi_{1i}) u(\xi_{1i}) w_i \) where \( outarray[p] = I[p], inarray[i] = u(\xi_{1i}), base[p*nq+i] = \phi_p(\xi_{1i}) \).

Parameters
basean array defining the local basis for the inner product usually passed from Basis->GetBdata() or Basis->GetDbdata()
inarrayphysical point array of function to be integrated \( u(\xi_1) \)
coll_checkflag to identify when a Basis->Collocation() call should be performed to see if this is a GLL_Lagrange basis with a collocation property. (should be set to 0 if taking the inner product with respect to the derivative of basis)
outarraythe values of the inner product with respect to each basis over region will be stored in the array outarray as output of the function

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 396 of file StdSegExp.cpp.

400{
401 boost::ignore_unused(coll_check);
402
403 int nquad = m_base[0]->GetNumPoints();
404 Array<OneD, NekDouble> tmp(nquad);
405 Array<OneD, const NekDouble> w = m_base[0]->GetW();
406
407 Vmath::Vmul(nquad, inarray, 1, w, 1, tmp, 1);
408
409 /* Comment below was a bug for collocated basis
410 if(coll_check&&m_base[0]->Collocation())
411 {
412 Vmath::Vcopy(nquad, tmp, 1, outarray, 1);
413 }
414 else
415 {
416 Blas::Dgemv('T',nquad,m_ncoeffs,1.0,base.get(),nquad,
417 &tmp[0],1,0.0,outarray.get(),1);
418 }*/
419
420 // Correct implementation
421 Blas::Dgemv('T', nquad, m_ncoeffs, 1.0, base.get(), nquad, &tmp[0], 1, 0.0,
422 outarray.get(), 1);
423}
std::vector< double > w(NPUPPER)

References Blas::Dgemv(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Vmath::Vmul(), and Nektar::UnitTests::w().

◆ v_IProductWRTBase() [2/2]

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

Inner product of inarray over region with respect to the expansion basis (this)->m_base[0] and return in outarray.

Wrapper call to StdSegExp::IProductWRTBase

Parameters
inarrayarray of function values evaluated at the physical collocation points
outarraythe values of the inner product with respect to each basis over region will be stored in the array outarray as output of the function

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 435 of file StdSegExp.cpp.

437{
438 v_IProductWRTBase(m_base[0]->GetBdata(), inarray, outarray, 1);
439}

References Nektar::StdRegions::StdExpansion::m_base, and v_IProductWRTBase().

Referenced by v_FwdTrans(), v_FwdTransBndConstrained(), v_HelmholtzMatrixOp(), v_IProductWRTBase(), v_IProductWRTDerivBase_SumFac(), and v_LaplacianMatrixOp().

◆ v_IProductWRTBase_SumFac()

void Nektar::StdRegions::StdSegExp::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 457 of file StdSegExp.cpp.

460{
461 int nquad = m_base[0]->GetNumPoints();
462 Array<OneD, NekDouble> tmp(nquad);
463 Array<OneD, const NekDouble> w = m_base[0]->GetW();
464 Array<OneD, const NekDouble> base = m_base[0]->GetBdata();
465
466 if (multiplybyweights)
467 {
468 Vmath::Vmul(nquad, inarray, 1, w, 1, tmp, 1);
469
470 Blas::Dgemv('T', nquad, m_ncoeffs, 1.0, base.get(), nquad, &tmp[0], 1,
471 0.0, outarray.get(), 1);
472 }
473 else
474 {
475 Blas::Dgemv('T', nquad, m_ncoeffs, 1.0, base.get(), nquad, &inarray[0],
476 1, 0.0, outarray.get(), 1);
477 }
478}

References Blas::Dgemv(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Vmath::Vmul(), and Nektar::UnitTests::w().

◆ v_IProductWRTDerivBase()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 441 of file StdSegExp.cpp.

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

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

◆ v_IProductWRTDerivBase_SumFac()

void Nektar::StdRegions::StdSegExp::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 448 of file StdSegExp.cpp.

451{
452 boost::ignore_unused(dir);
453 ASSERTL1(dir == 0, "input dir is out of range");
454 StdSegExp::v_IProductWRTBase(m_base[0]->GetDbdata(), inarray, outarray, 1);
455}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249

References ASSERTL1, Nektar::StdRegions::StdExpansion::m_base, and v_IProductWRTBase().

◆ v_IsBoundaryInteriorExpansion()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 81 of file StdSegExp.cpp.

82{
83
84 bool returnval = false;
85
87 {
88 returnval = true;
89 }
90
92 {
93 returnval = true;
94 }
95
96 return returnval;
97}

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

◆ v_LaplacianMatrixOp()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 512 of file StdSegExp.cpp.

515{
516 boost::ignore_unused(mkey);
517
518 int nquad = m_base[0]->GetNumPoints();
519
520 Array<OneD, NekDouble> physValues(nquad);
521 Array<OneD, NekDouble> dPhysValuesdx(nquad);
522
523 v_BwdTrans(inarray, physValues);
524
525 // Laplacian matrix operation
526 v_PhysDeriv(physValues, dPhysValuesdx);
527 v_IProductWRTBase(m_base[0]->GetDbdata(), dPhysValuesdx, outarray, 1);
528}

References Nektar::StdRegions::StdExpansion::m_base, v_BwdTrans(), v_IProductWRTBase(), and v_PhysDeriv().

◆ v_LocCollapsedToLocCoord()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 489 of file StdSegExp.cpp.

491{
492 xi[0] = eta[0];
493}

◆ v_LocCoordToLocCollapsed()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 483 of file StdSegExp.cpp.

485{
486 eta[0] = xi[0];
487}

◆ v_MultiplyByStdQuadratureMetric()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 633 of file StdSegExp.cpp.

636{
637 int nquad0 = m_base[0]->GetNumPoints();
638
639 const Array<OneD, const NekDouble> &w0 = m_base[0]->GetW();
640
641 Vmath::Vmul(nquad0, inarray.get(), 1, w0.get(), 1, outarray.get(), 1);
642}

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

◆ v_NumBndryCoeffs()

int Nektar::StdRegions::StdSegExp::v_NumBndryCoeffs ( ) const
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 684 of file StdSegExp.cpp.

685{
686 return 2;
687}

◆ v_NumDGBndryCoeffs()

int Nektar::StdRegions::StdSegExp::v_NumDGBndryCoeffs ( ) const
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 689 of file StdSegExp.cpp.

690{
691 return 2;
692}

◆ v_PhysDeriv() [1/2]

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

Evaluate the derivative \( d/d{\xi_1} \) at the physical quadrature points given by inarray and return in outarray.

This is a wrapper around StdExpansion1D::Tensor_Deriv

Parameters
inarrayarray of a function evaluated at the quadrature points
outarraythe resulting array of the derivative \( du/d_{\xi_1}|_{\xi_{1i}} \) will be stored in the array outarra

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 140 of file StdSegExp.cpp.

144{
145 boost::ignore_unused(out_d1, out_d2);
146 PhysTensorDeriv(inarray, out_d0);
147}
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Evaluate the derivative at the physical quadrature points given by inarray and return in outarray.

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

Referenced by v_HelmholtzMatrixOp(), and v_LaplacianMatrixOp().

◆ v_PhysDeriv() [2/2]

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

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

See also
StdRegions::StdExpansion::PhysDeriv

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 149 of file StdSegExp.cpp.

152{
153 boost::ignore_unused(dir);
154 ASSERTL1(dir == 0, "input dir is out of range");
155 PhysTensorDeriv(inarray, outarray);
156 // PhysDeriv(inarray, outarray);
157}

References ASSERTL1, and Nektar::StdRegions::StdExpansion1D::PhysTensorDeriv().

◆ v_PhysEvaluate() [1/2]

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

Reimplemented from Nektar::StdRegions::StdExpansion1D.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 150 of file StdSegExp.h.

154 {
155 return StdExpansion1D::BaryTensorDeriv(coord, inarray,
156 firstOrderDerivs);
157 }
NekDouble BaryTensorDeriv(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)

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

◆ v_PhysEvaluate() [2/2]

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

Reimplemented from Nektar::StdRegions::StdExpansion1D.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 159 of file StdSegExp.h.

164 {
165 return StdExpansion1D::BaryTensorDeriv(coord, inarray, firstOrderDerivs,
166 secondOrderDerivs);
167 }

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

◆ v_PhysEvaluateBasis()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 506 of file StdSegExp.cpp.

508{
509 return StdExpansion::BaryEvaluateBasis<0>(coords[0], mode);
510}

◆ v_ReduceOrderCoeffs()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 220 of file StdSegExp.cpp.

223{
224 int n_coeffs = inarray.size();
225
226 Array<OneD, NekDouble> coeff(n_coeffs);
227 Array<OneD, NekDouble> coeff_tmp(n_coeffs, 0.0);
228 Array<OneD, NekDouble> tmp;
229 Array<OneD, NekDouble> tmp2;
230
231 int nmodes0 = m_base[0]->GetNumModes();
232
233 Vmath::Vcopy(n_coeffs, inarray, 1, coeff_tmp, 1);
234
235 const LibUtilities::PointsKey Pkey0(nmodes0,
237
238 LibUtilities::BasisKey b0(m_base[0]->GetBasisType(), nmodes0, Pkey0);
239
240 LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A, nmodes0, Pkey0);
241
242 LibUtilities::InterpCoeff1D(b0, coeff_tmp, bortho0, coeff);
243
244 Vmath::Zero(n_coeffs, coeff_tmp, 1);
245
246 Vmath::Vcopy(numMin, tmp = coeff, 1, tmp2 = coeff_tmp, 1);
247
248 LibUtilities::InterpCoeff1D(bortho0, coeff_tmp, b0, outarray);
249}
void InterpCoeff1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
Definition: InterpCoeff.cpp:44
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:53
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487

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

◆ v_StdPhysDeriv() [1/2]

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 159 of file StdSegExp.cpp.

163{
164 boost::ignore_unused(out_d1, out_d2);
165 PhysTensorDeriv(inarray, out_d0);
166 // PhysDeriv(inarray, out_d0);
167}

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

◆ v_StdPhysDeriv() [2/2]

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 169 of file StdSegExp.cpp.

172{
173 boost::ignore_unused(dir);
174 ASSERTL1(dir == 0, "input dir is out of range");
175 PhysTensorDeriv(inarray, outarray);
176 // PhysDeriv(inarray, outarray);
177}

References ASSERTL1, and Nektar::StdRegions::StdExpansion1D::PhysTensorDeriv().

◆ v_SVVLaplacianFilter()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 552 of file StdSegExp.cpp.

554{
555 // Generate an orthogonal expansion
556 int nq = m_base[0]->GetNumPoints();
557 int nmodes = m_base[0]->GetNumModes();
558 // Declare orthogonal basis.
559 LibUtilities::PointsKey pKey(nq, m_base[0]->GetPointsType());
560
561 LibUtilities::BasisKey B(LibUtilities::eOrtho_A, nmodes, pKey);
562 StdSegExp OrthoExp(B);
563
564 // SVV parameters loaded from the .xml case file
565 NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDiffCoeff);
566 int cutoff = (int)(mkey.GetConstFactor(eFactorSVVCutoffRatio)) * nmodes;
567
568 Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
569
570 // project onto modal space.
571 OrthoExp.FwdTrans(array, orthocoeffs);
572
573 //
574 for (int j = 0; j < nmodes; ++j)
575 {
576 if (j >= cutoff) // to filter out only the "high-modes"
577 {
578 orthocoeffs[j] *=
579 (SvvDiffCoeff *
580 exp(-(j - nmodes) * (j - nmodes) /
581 ((NekDouble)((j - cutoff + 1) * (j - cutoff + 1)))));
582 }
583 else
584 {
585 orthocoeffs[j] *= 0.0;
586 }
587 }
588
589 // backward transform to physical space
590 OrthoExp.BwdTrans(orthocoeffs, array);
591}

References Nektar::StdRegions::StdExpansion::BwdTrans(), Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::eFactorSVVDiffCoeff, Nektar::LibUtilities::eOrtho_A, Nektar::StdRegions::StdExpansion::FwdTrans(), Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), Nektar::StdRegions::StdExpansion::GetPointsType(), and Nektar::StdRegions::StdExpansion::m_base.