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 (const LibUtilities::BasisKey &Ba)
 Constructor using BasisKey class for quadrature points and order definition. More...
 
 StdSegExp (const StdSegExp &T)=default
 
 ~StdSegExp () override=default
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion1D
 StdExpansion1D (int numcoeffs, const LibUtilities::BasisKey &Ba)
 
 StdExpansion1D ()=default
 
 StdExpansion1D (const StdExpansion1D &T)=default
 
 ~StdExpansion1D () override=default
 
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

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_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_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_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
 
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
 
NekDouble v_PhysEvaluateBasis (const Array< OneD, const NekDouble > &coords, int mode) final
 
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
 
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
 
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_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_FillMode (const int mode, Array< OneD, NekDouble > &outarray) override
 
void v_GetCoords (Array< OneD, NekDouble > &coords_0, Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2) 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
 
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 override
 
int v_NumDGBndryCoeffs () const override
 
bool v_IsBoundaryInteriorExpansion () const override
 
int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset) override
 
LibUtilities::ShapeType v_DetShapeType () const override
 Return Shape of region, using ShapeType enum list. i.e. Segment. More...
 
DNekMatSharedPtr v_GenMatrix (const StdMatrixKey &mkey) override
 
DNekMatSharedPtr v_CreateStdMatrix (const StdMatrixKey &mkey) override
 
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
 
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
NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals) override
 
NekDouble v_PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
 
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 44 of file StdSegExp.h.

Constructor & Destructor Documentation

◆ StdSegExp() [1/2]

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 48 of file StdSegExp.cpp.

49 : StdExpansion(Ba.GetNumModes(), 1, Ba),
50 StdExpansion1D(Ba.GetNumModes(), Ba)
51{
52}
StdExpansion()
Default Constructor.

◆ StdSegExp() [2/2]

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

◆ ~StdSegExp()

Nektar::StdRegions::StdSegExp::~StdSegExp ( )
overridedefault

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 177 of file StdSegExp.cpp.

179{
180 int nquad = m_base[0]->GetNumPoints();
181
182 if (m_base[0]->Collocation())
183 {
184 Vmath::Vcopy(nquad, inarray, 1, outarray, 1);
185 }
186 else
187 {
188 Blas::Dgemv('N', nquad, m_base[0]->GetNumModes(), 1.0,
189 (m_base[0]->GetBdata()).get(), nquad, &inarray[0], 1, 0.0,
190 &outarray[0], 1);
191 }
192}
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:211
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

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 339 of file StdSegExp.cpp.

341{
342 v_BwdTrans(inarray, outarray);
343}
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:177

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 659 of file StdSegExp.cpp.

661{
662 int nmodes = nummodes[modes_offset];
663 modes_offset += 1;
664
665 return nmodes;
666}

◆ 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 735 of file StdSegExp.cpp.

736{
737 return v_GenMatrix(mkey);
738}
DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override
Definition: StdSegExp.cpp:672

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 57 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 562 of file StdSegExp.cpp.

566{
567 // Generate an orthogonal expansion
568 int nq = m_base[0]->GetNumPoints();
569 int nmodes = m_base[0]->GetNumModes();
570 int P = nmodes - 1;
571 // Declare orthogonal basis.
572 LibUtilities::PointsKey pKey(nq, m_base[0]->GetPointsType());
573
574 LibUtilities::BasisKey B(LibUtilities::eOrtho_A, nmodes, pKey);
575 StdSegExp OrthoExp(B);
576
577 // Cutoff
578 int Pcut = cutoff * P;
579
580 // Project onto orthogonal space.
581 Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
582 OrthoExp.FwdTrans(array, orthocoeffs);
583
584 //
585 NekDouble fac;
586 for (int j = 0; j < nmodes; ++j)
587 {
588 // to filter out only the "high-modes"
589 if (j > Pcut)
590 {
591 fac = (NekDouble)(j - Pcut) / ((NekDouble)(P - Pcut));
592 fac = pow(fac, exponent);
593 orthocoeffs[j] *= exp(-alpha * fac);
594 }
595 }
596
597 // backward transform to physical space
598 OrthoExp.BwdTrans(orthocoeffs, array);
599}
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:205
StdSegExp(const LibUtilities::BasisKey &Ba)
Constructor using BasisKey class for quadrature points and order definition.
Definition: StdSegExp.cpp:48
@ P
Monomial polynomials .
Definition: BasisType.h:62
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:42
double NekDouble

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

◆ v_FillMode()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 466 of file StdSegExp.cpp.

467{
468 int nquad = m_base[0]->GetNumPoints();
469 const NekDouble *base = m_base[0]->GetBdata().get();
470
471 ASSERTL2(mode <= m_ncoeffs,
472 "calling argument mode is larger than total expansion order");
473
474 Vmath::Vcopy(nquad, (NekDouble *)base + mode * nquad, 1, &outarray[0], 1);
475}
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:265

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 243 of file StdSegExp.cpp.

245{
246 if (m_base[0]->Collocation())
247 {
248 Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
249 }
250 else
251 {
252 v_IProductWRTBase(inarray, outarray);
253
254 // get Mass matrix inverse
255 StdMatrixKey masskey(eInvMass, v_DetShapeType(), *this);
256 DNekMatSharedPtr matsys = GetStdMatrix(masskey);
257
258 NekVector<NekDouble> in(m_ncoeffs, outarray, eCopy);
259 NekVector<NekDouble> out(m_ncoeffs, outarray, eWrapper);
260
261 out = (*matsys) * in;
262 }
263}
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:603
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:407
LibUtilities::ShapeType v_DetShapeType() const override
Return Shape of region, using ShapeType enum list. i.e. Segment.
Definition: StdSegExp.cpp:57
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 265 of file StdSegExp.cpp.

268{
269 if (m_base[0]->Collocation())
270 {
271 Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
272 }
273 else
274 {
275 int nInteriorDofs = m_ncoeffs - 2;
276 int offset = 0;
277
278 switch (m_base[0]->GetBasisType())
279 {
281 {
282 offset = 1;
283 }
284 break;
286 {
287 nInteriorDofs = m_ncoeffs;
288 offset = 0;
289 }
290 break;
293 {
294 offset = 2;
295 }
296 break;
297 default:
298 ASSERTL0(false, "This type of FwdTrans is not defined for this "
299 "expansion type");
300 }
301
302 fill(outarray.get(), outarray.get() + m_ncoeffs, 0.0);
303
305 {
306 outarray[GetVertexMap(0)] = inarray[0];
307 outarray[GetVertexMap(1)] = inarray[m_base[0]->GetNumPoints() - 1];
308
309 if (m_ncoeffs > 2)
310 {
311 // ideally, we would like to have tmp0 to be replaced by
312 // outarray (currently MassMatrixOp does not allow aliasing)
313 Array<OneD, NekDouble> tmp0(m_ncoeffs);
314 Array<OneD, NekDouble> tmp1(m_ncoeffs);
315
316 StdMatrixKey masskey(eMass, v_DetShapeType(), *this);
317 MassMatrixOp(outarray, tmp0, masskey);
318 v_IProductWRTBase(inarray, tmp1);
319
320 Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
321
322 // get Mass matrix inverse (only of interior DOF)
323 DNekMatSharedPtr matsys =
324 (m_stdStaticCondMatrixManager[masskey])->GetBlock(1, 1);
325
326 Blas::Dgemv('N', nInteriorDofs, nInteriorDofs, 1.0,
327 &(matsys->GetPtr())[0], nInteriorDofs,
328 tmp1.get() + offset, 1, 0.0,
329 outarray.get() + offset, 1);
330 }
331 }
332 else
333 {
334 StdSegExp::v_FwdTrans(inarray, outarray);
335 }
336 }
337}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:156
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:752
int GetVertexMap(const int localVertexId, bool useCoeffPacking=false)
Definition: StdExpansion.h:679
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:243
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:49
@ eGauss_Lagrange
Lagrange Polynomials using the Gauss points.
Definition: BasisType.h:57
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:56
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48
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 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 672 of file StdSegExp.cpp.

673{
675 MatrixType mattype;
676
677 switch (mattype = mkey.GetMatrixType())
678 {
680 {
681 int nq = m_base[0]->GetNumPoints();
682
683 // take definition from key
684 if (mkey.ConstFactorExists(eFactorConst))
685 {
686 nq = (int)mkey.GetConstFactor(eFactorConst);
687 }
688
690 Array<OneD, NekDouble> coords(1);
693
694 for (int i = 0; i < neq; ++i)
695 {
696 coords[0] = -1.0 + 2 * i / (NekDouble)(neq - 1);
697 I = m_base[0]->GetI(coords);
698 Vmath::Vcopy(nq, I->GetRawPtr(), 1, Mat->GetRawPtr() + i, neq);
699 }
700 }
701 break;
702 case eFwdTrans:
703 {
704 Mat =
706 StdMatrixKey iprodkey(eIProductWRTBase, v_DetShapeType(), *this);
707 DNekMat &Iprod = *GetStdMatrix(iprodkey);
708 StdMatrixKey imasskey(eInvMass, v_DetShapeType(), *this);
709 DNekMat &Imass = *GetStdMatrix(imasskey);
710
711 (*Mat) = Imass * Iprod;
712 }
713 break;
714 default:
715 {
717
718 if (mattype == eMass)
719 {
720 // For Fourier basis set the imaginary component
721 // of mean mode to have a unit diagonal component
722 // in mass matrix
724 {
725 (*Mat)(1, 1) = 1.0;
726 }
727 }
728 }
729 break;
730 }
731
732 return Mat;
733}
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:55
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 744 of file StdSegExp.cpp.

745{
746 if (outarray.size() != NumBndryCoeffs())
747 {
748 outarray = Array<OneD, unsigned int>(NumBndryCoeffs());
749 }
750 const LibUtilities::BasisType Btype = GetBasisType(0);
751 int nummodes = m_base[0]->GetNumModes();
752
753 outarray[0] = 0;
754
755 switch (Btype)
756 {
761 outarray[1] = nummodes - 1;
762 break;
765 outarray[1] = 1;
766 break;
767 default:
768 ASSERTL0(0, "Mapping array is not defined for this expansion");
769 break;
770 }
771}
@ eChebyshev
Chebyshev Polynomials.
Definition: BasisType.h:61

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 613 of file StdSegExp.cpp.

616{
617 Blas::Dcopy(GetNumPoints(0), (m_base[0]->GetZ()).get(), 1, &coords_0[0], 1);
618}
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:218
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

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 889 of file StdSegExp.cpp.

894{
895 // parameters for higher dimnesion traces
896 if (maparray.size() != 1)
897 {
898 maparray = Array<OneD, unsigned int>(1);
899 }
900
901 maparray[0] = 0;
902
903 if (signarray.size() != 1)
904 {
905 signarray = Array<OneD, int>(1, 1);
906 }
907 else
908 {
909 signarray[0] = 1;
910 }
911}

◆ v_GetInteriorMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 773 of file StdSegExp.cpp.

774{
775 int i;
776 if (outarray.size() != GetNcoeffs() - NumBndryCoeffs())
777 {
778 outarray = Array<OneD, unsigned int>(GetNcoeffs() - NumBndryCoeffs());
779 }
780 const LibUtilities::BasisType Btype = GetBasisType(0);
781
782 switch (Btype)
783 {
788 for (i = 0; i < GetNcoeffs() - 2; i++)
789 {
790 outarray[i] = i + 1;
791 }
792 break;
795 for (i = 0; i < GetNcoeffs() - 2; i++)
796 {
797 outarray[i] = i + 2;
798 }
799 break;
800 default:
801 ASSERTL0(0, "Mapping array is not defined for this expansion");
802 break;
803 }
804}
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:124

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
finalprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 629 of file StdSegExp.cpp.

630{
631 return 2;
632}

◆ v_GetNverts()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 624 of file StdSegExp.cpp.

625{
626 return 2;
627}

◆ v_GetSimplexEquiSpacedConnectivity()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 823 of file StdSegExp.cpp.

825{
826 int np = m_base[0]->GetNumPoints();
827
828 conn = Array<OneD, int>(2 * (np - 1));
829 int cnt = 0;
830 for (int i = 0; i < np - 1; ++i)
831 {
832 conn[cnt++] = i;
833 conn[cnt++] = i + 1;
834 }
835}

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 841 of file StdSegExp.cpp.

843{
844 int order0 = m_base[0]->GetNumModes();
845
846 ASSERTL0(traceid < 2, "eid must be between 0 and 1");
847
848 if (maparray.size() != 1)
849 {
850 maparray = Array<OneD, unsigned int>(1);
851 }
852
853 const LibUtilities::BasisType bType = GetBasisType(0);
854
855 if (bType == LibUtilities::eModified_A)
856 {
857 maparray[0] = (traceid == 0) ? 0 : 1;
858 }
859 else if (bType == LibUtilities::eGLL_Lagrange ||
861 {
862 maparray[0] = (traceid == 0) ? 0 : order0 - 1;
863 }
864 else
865 {
866 ASSERTL0(false, "Unknown Basis");
867 }
868}

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
finalprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 639 of file StdSegExp.cpp.

640{
641 return 0;
642}

◆ v_GetTraceNcoeffs()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 634 of file StdSegExp.cpp.

635{
636 return 1;
637}

◆ v_GetTraceNumPoints()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 644 of file StdSegExp.cpp.

645{
646 return 1;
647}

◆ 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 870 of file StdSegExp.cpp.

876{
877 v_GetTraceCoeffMap(tid, maparray);
878
879 if (signarray.size() != 1)
880 {
881 signarray = Array<OneD, int>(1, 1);
882 }
883 else
884 {
885 signarray[0] = 1;
886 }
887}
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:841

References v_GetTraceCoeffMap().

◆ v_GetVertexMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 806 of file StdSegExp.cpp.

808{
809 ASSERTL0((localVertexId == 0) || (localVertexId == 1),
810 "local vertex id"
811 "must be between 0 or 1");
812
813 int localDOF = localVertexId;
814
816 (localVertexId == 1))
817 {
818 localDOF = m_base[0]->GetNumModes() - 1;
819 }
820 return localDOF;
821}

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 499 of file StdSegExp.cpp.

502{
503 int nquad = m_base[0]->GetNumPoints();
504
505 Array<OneD, NekDouble> physValues(nquad);
506 Array<OneD, NekDouble> dPhysValuesdx(nquad);
507 Array<OneD, NekDouble> wsp(m_ncoeffs);
508
509 v_BwdTrans(inarray, physValues);
510
511 // mass matrix operation
512 v_IProductWRTBase((m_base[0]->GetBdata()), physValues, wsp, 1);
513
514 // Laplacian matrix operation
515 v_PhysDeriv(physValues, dPhysValuesdx);
516 v_IProductWRTBase(m_base[0]->GetDbdata(), dPhysValuesdx, outarray, 1);
517 Blas::Daxpy(m_ncoeffs, mkey.GetConstFactor(eFactorLambda), wsp.get(), 1,
518 outarray.get(), 1);
519}
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:121
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:135

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 92 of file StdSegExp.cpp.

93{
94 NekDouble Int = 0.0;
95 int nquad0 = m_base[0]->GetNumPoints();
96 Array<OneD, NekDouble> tmp(nquad0);
97 Array<OneD, const NekDouble> z = m_base[0]->GetZ();
98 Array<OneD, const NekDouble> w0 = m_base[0]->GetW();
99
100 // multiply by integration constants
101 Vmath::Vmul(nquad0, inarray, 1, w0, 1, tmp, 1);
102
103 Int = Vmath::Vsum(nquad0, tmp, 1);
104
105 return Int;
106}
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.hpp:72
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.hpp:608

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 370 of file StdSegExp.cpp.

374{
375 int nquad = m_base[0]->GetNumPoints();
376 Array<OneD, NekDouble> tmp(nquad);
377 Array<OneD, const NekDouble> w = m_base[0]->GetW();
378
379 Vmath::Vmul(nquad, inarray, 1, w, 1, tmp, 1);
380
381 /* Comment below was a bug for collocated basis
382 if(coll_check&&m_base[0]->Collocation())
383 {
384 Vmath::Vcopy(nquad, tmp, 1, outarray, 1);
385 }
386 else
387 {
388 Blas::Dgemv('T',nquad,m_ncoeffs,1.0,base.get(),nquad,
389 &tmp[0],1,0.0,outarray.get(),1);
390 }*/
391
392 // Correct implementation
393 Blas::Dgemv('T', nquad, m_ncoeffs, 1.0, base.get(), nquad, &tmp[0], 1, 0.0,
394 outarray.get(), 1);
395}
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 407 of file StdSegExp.cpp.

409{
410 v_IProductWRTBase(m_base[0]->GetBdata(), inarray, outarray, 1);
411}

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 428 of file StdSegExp.cpp.

431{
432 int nquad = m_base[0]->GetNumPoints();
433 Array<OneD, NekDouble> tmp(nquad);
434 Array<OneD, const NekDouble> w = m_base[0]->GetW();
435 Array<OneD, const NekDouble> base = m_base[0]->GetBdata();
436
437 if (multiplybyweights)
438 {
439 Vmath::Vmul(nquad, inarray, 1, w, 1, tmp, 1);
440
441 Blas::Dgemv('T', nquad, m_ncoeffs, 1.0, base.get(), nquad, &tmp[0], 1,
442 0.0, outarray.get(), 1);
443 }
444 else
445 {
446 Blas::Dgemv('T', nquad, m_ncoeffs, 1.0, base.get(), nquad, &inarray[0],
447 1, 0.0, outarray.get(), 1);
448 }
449}

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 413 of file StdSegExp.cpp.

416{
417 StdSegExp::IProductWRTDerivBase_SumFac(dir, inarray, outarray);
418}
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 420 of file StdSegExp.cpp.

423{
424 ASSERTL1(dir == 0, "input dir is out of range");
425 StdSegExp::v_IProductWRTBase(m_base[0]->GetDbdata(), inarray, outarray, 1);
426}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242

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 62 of file StdSegExp.cpp.

63{
64
65 bool returnval = false;
66
68 {
69 returnval = true;
70 }
71
73 {
74 returnval = true;
75 }
76
77 return returnval;
78}

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 483 of file StdSegExp.cpp.

486{
487 int nquad = m_base[0]->GetNumPoints();
488
489 Array<OneD, NekDouble> physValues(nquad);
490 Array<OneD, NekDouble> dPhysValuesdx(nquad);
491
492 v_BwdTrans(inarray, physValues);
493
494 // Laplacian matrix operation
495 v_PhysDeriv(physValues, dPhysValuesdx);
496 v_IProductWRTBase(m_base[0]->GetDbdata(), dPhysValuesdx, outarray, 1);
497}

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 460 of file StdSegExp.cpp.

462{
463 xi[0] = eta[0];
464}

◆ 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 454 of file StdSegExp.cpp.

456{
457 eta[0] = xi[0];
458}

◆ 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 602 of file StdSegExp.cpp.

605{
606 int nquad0 = m_base[0]->GetNumPoints();
607
608 const Array<OneD, const NekDouble> &w0 = m_base[0]->GetW();
609
610 Vmath::Vmul(nquad0, inarray.get(), 1, w0.get(), 1, outarray.get(), 1);
611}

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 649 of file StdSegExp.cpp.

650{
651 return 2;
652}

◆ v_NumDGBndryCoeffs()

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

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 654 of file StdSegExp.cpp.

655{
656 return 2;
657}

◆ 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 121 of file StdSegExp.cpp.

125{
126 PhysTensorDeriv(inarray, out_d0);
127}
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 129 of file StdSegExp.cpp.

132{
133 ASSERTL1(dir == 0, "input dir is out of range");
134 PhysTensorDeriv(inarray, outarray);
135 // PhysDeriv(inarray, outarray);
136}

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 126 of file StdSegExp.h.

130 {
131 return StdExpansion1D::BaryTensorDeriv(coord, inarray,
132 firstOrderDerivs);
133 }
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 135 of file StdSegExp.h.

139 {
140 return StdExpansion1D::BaryTensorDeriv(coord, inarray, firstOrderDerivs,
141 secondOrderDerivs);
142 }

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

◆ v_PhysEvaluateBasis()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 477 of file StdSegExp.cpp.

479{
480 return StdExpansion::BaryEvaluateBasis<0>(coords[0], mode);
481}

◆ 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 194 of file StdSegExp.cpp.

197{
198 int n_coeffs = inarray.size();
199
200 Array<OneD, NekDouble> coeff(n_coeffs);
201 Array<OneD, NekDouble> coeff_tmp(n_coeffs, 0.0);
202 Array<OneD, NekDouble> tmp;
203 Array<OneD, NekDouble> tmp2;
204
205 int nmodes0 = m_base[0]->GetNumModes();
206
207 Vmath::Vcopy(n_coeffs, inarray, 1, coeff_tmp, 1);
208
209 const LibUtilities::PointsKey Pkey0(nmodes0,
211
212 LibUtilities::BasisKey b0(m_base[0]->GetBasisType(), nmodes0, Pkey0);
213
214 LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A, nmodes0, Pkey0);
215
216 LibUtilities::InterpCoeff1D(b0, coeff_tmp, bortho0, coeff);
217
218 Vmath::Zero(n_coeffs, coeff_tmp, 1);
219
220 Vmath::Vcopy(numMin, tmp = coeff, 1, tmp2 = coeff_tmp, 1);
221
222 LibUtilities::InterpCoeff1D(bortho0, coeff_tmp, b0, outarray);
223}
void InterpCoeff1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
Definition: InterpCoeff.cpp:42
@ 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::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 138 of file StdSegExp.cpp.

142{
143 PhysTensorDeriv(inarray, out_d0);
144}

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 146 of file StdSegExp.cpp.

149{
150 ASSERTL1(dir == 0, "input dir is out of range");
151 PhysTensorDeriv(inarray, outarray);
152}

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 521 of file StdSegExp.cpp.

523{
524 // Generate an orthogonal expansion
525 int nq = m_base[0]->GetNumPoints();
526 int nmodes = m_base[0]->GetNumModes();
527 // Declare orthogonal basis.
528 LibUtilities::PointsKey pKey(nq, m_base[0]->GetPointsType());
529
530 LibUtilities::BasisKey B(LibUtilities::eOrtho_A, nmodes, pKey);
531 StdSegExp OrthoExp(B);
532
533 // SVV parameters loaded from the .xml case file
534 NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDiffCoeff);
535 int cutoff = (int)(mkey.GetConstFactor(eFactorSVVCutoffRatio)) * nmodes;
536
537 Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
538
539 // project onto modal space.
540 OrthoExp.FwdTrans(array, orthocoeffs);
541
542 //
543 for (int j = 0; j < nmodes; ++j)
544 {
545 if (j >= cutoff) // to filter out only the "high-modes"
546 {
547 orthocoeffs[j] *=
548 (SvvDiffCoeff *
549 exp(-(j - nmodes) * (j - nmodes) /
550 ((NekDouble)((j - cutoff + 1) * (j - cutoff + 1)))));
551 }
552 else
553 {
554 orthocoeffs[j] *= 0.0;
555 }
556 }
557
558 // backward transform to physical space
559 OrthoExp.BwdTrans(orthocoeffs, array);
560}

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.