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

Protected Member Functions

NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray) override
 Integrate the physical point list inarray over region and return the value.
 
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.
 
void v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Calculate the derivative of the physical points in a given direction.
 
void v_StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray) 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.
 
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.
 
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.
 
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.
 
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_PhysEvalFirstDeriv (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
 
NekDouble v_PhysEvalFirstSecondDeriv (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_LaplacianMatrixOp (const int k1, const int k2, 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.
 
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.
 
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
 
void v_PhysInterp (std::shared_ptr< StdExpansion > fromExp, const Array< OneD, const NekDouble > &fromData, Array< OneD, NekDouble > &toData) override
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion
DNekMatSharedPtr CreateStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr CreateStdStaticCondMatrix (const StdMatrixKey &mkey)
 Create the static condensation of a matrix when using a boundary interior decomposition.
 
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase_SumFac (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GeneralMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
 
void LaplacianMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp_MatFree (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void HelmholtzMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
virtual void v_GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals, NekDouble &deriv, NekDouble &deriv2)
 This function performs the barycentric interpolation of the polynomial stored in coord at a point physvals using barycentric interpolation weights in direction.
 
template<int DIR>
NekDouble BaryEvaluateBasis (const NekDouble &coord, const int &mode)
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals)
 Helper function to pass an unused value by reference into BaryEvaluate.
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals, NekDouble &deriv)
 

Additional Inherited Members

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

Detailed Description

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

171{
172 int nquad = m_base[0]->GetNumPoints();
173
174 if (m_base[0]->Collocation())
175 {
176 Vmath::Vcopy(nquad, inarray, 1, outarray, 1);
177 }
178 else
179 {
180 Blas::Dgemv('N', nquad, m_base[0]->GetNumModes(), 1.0,
181 (m_base[0]->GetBdata()).data(), nquad, &inarray[0], 1, 0.0,
182 &outarray[0], 1);
183 }
184}
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 331 of file StdSegExp.cpp.

333{
334 v_BwdTrans(inarray, outarray);
335}
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...

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

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

◆ v_CreateStdMatrix()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 753 of file StdSegExp.cpp.

754{
755 return v_GenMatrix(mkey);
756}
DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override

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

583{
584 // Generate an orthogonal expansion
585 int nq = m_base[0]->GetNumPoints();
586 int nmodes = m_base[0]->GetNumModes();
587 int P = nmodes - 1;
588 // Declare orthogonal basis.
589 LibUtilities::PointsKey pKey(nq, m_base[0]->GetPointsType());
590
591 LibUtilities::BasisKey B(LibUtilities::eOrtho_A, nmodes, pKey);
592 StdSegExp OrthoExp(B);
593
594 // Cutoff
595 int Pcut = cutoff * P;
596
597 // Project onto orthogonal space.
598 Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
599 OrthoExp.FwdTrans(array, orthocoeffs);
600
601 //
602 NekDouble fac;
603 for (int j = 0; j < nmodes; ++j)
604 {
605 // to filter out only the "high-modes"
606 if (j > Pcut)
607 {
608 fac = (NekDouble)(j - Pcut) / ((NekDouble)(P - Pcut));
609 fac = pow(fac, exponent);
610 orthocoeffs[j] *= exp(-alpha * fac);
611 }
612 }
613
614 // backward transform to physical space
615 OrthoExp.BwdTrans(orthocoeffs, array);
616}
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
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

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

459{
460 int nquad = m_base[0]->GetNumPoints();
461 const NekDouble *base = m_base[0]->GetBdata().data();
462
463 ASSERTL2(mode <= m_ncoeffs,
464 "calling argument mode is larger than total expansion order");
465
466 Vmath::Vcopy(nquad, (NekDouble *)base + mode * nquad, 1, &outarray[0], 1);
467}
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
const Array< OneD, NekDouble > base(3, 0.0)

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.

Definition at line 235 of file StdSegExp.cpp.

237{
238 if (m_base[0]->Collocation())
239 {
240 Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
241 }
242 else
243 {
244 v_IProductWRTBase(inarray, outarray);
245
246 // get Mass matrix inverse
247 StdMatrixKey masskey(eInvMass, v_DetShapeType(), *this);
248 DNekMatSharedPtr matsys = GetStdMatrix(masskey);
249
250 NekVector<NekDouble> in(m_ncoeffs, outarray, eCopy);
251 NekVector<NekDouble> out(m_ncoeffs, outarray, eWrapper);
252
253 out = (*matsys) * in;
254 }
255}
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
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...
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

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.

Definition at line 257 of file StdSegExp.cpp.

260{
261 if (m_base[0]->Collocation())
262 {
263 Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
264 }
265 else
266 {
267 int nInteriorDofs = m_ncoeffs - 2;
268 int offset = 0;
269
270 switch (m_base[0]->GetBasisType())
271 {
273 {
274 offset = 1;
275 }
276 break;
278 {
279 nInteriorDofs = m_ncoeffs;
280 offset = 0;
281 }
282 break;
285 {
286 offset = 2;
287 }
288 break;
289 default:
290 ASSERTL0(false, "This type of FwdTrans is not defined for this "
291 "expansion type");
292 }
293
294 fill(outarray.data(), outarray.data() + m_ncoeffs, 0.0);
295
297 {
298 outarray[GetVertexMap(0)] = inarray[0];
299 outarray[GetVertexMap(1)] = inarray[m_base[0]->GetNumPoints() - 1];
300
301 if (m_ncoeffs > 2)
302 {
303 // ideally, we would like to have tmp0 to be replaced by
304 // outarray (currently MassMatrixOp does not allow aliasing)
305 Array<OneD, NekDouble> tmp0(m_ncoeffs);
306 Array<OneD, NekDouble> tmp1(m_ncoeffs);
307
308 StdMatrixKey masskey(eMass, v_DetShapeType(), *this);
309 MassMatrixOp(outarray, tmp0, masskey);
310 v_IProductWRTBase(inarray, tmp1);
311
312 Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
313
314 // get Mass matrix inverse (only of interior DOF)
315 DNekMatSharedPtr matsys =
316 (m_stdStaticCondMatrixManager[masskey])->GetBlock(1, 1);
317
318 Blas::Dgemv('N', nInteriorDofs, nInteriorDofs, 1.0,
319 &(matsys->GetPtr())[0], nInteriorDofs,
320 tmp1.data() + offset, 1, 0.0,
321 outarray.data() + offset, 1);
322 }
323 }
324 else
325 {
326 StdSegExp::v_FwdTrans(inarray, outarray);
327 }
328 }
329}
#define ASSERTL0(condition, msg)
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
int GetVertexMap(const int localVertexId, bool useCoeffPacking=false)
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...
@ 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.

Definition at line 690 of file StdSegExp.cpp.

691{
693 MatrixType mattype;
694
695 switch (mattype = mkey.GetMatrixType())
696 {
698 {
699 int nq = m_base[0]->GetNumPoints();
700
701 // take definition from key
702 if (mkey.ConstFactorExists(eFactorConst))
703 {
704 nq = (int)mkey.GetConstFactor(eFactorConst);
705 }
706
708 Array<OneD, NekDouble> coords(1);
711
712 for (int i = 0; i < neq; ++i)
713 {
714 coords[0] = -1.0 + 2 * i / (NekDouble)(neq - 1);
715 I = m_base[0]->GetI(coords);
716 Vmath::Vcopy(nq, I->GetRawPtr(), 1, Mat->GetRawPtr() + i, neq);
717 }
718 }
719 break;
720 case eFwdTrans:
721 {
722 Mat =
724 StdMatrixKey iprodkey(eIProductWRTBase, v_DetShapeType(), *this);
725 DNekMat &Iprod = *GetStdMatrix(iprodkey);
726 StdMatrixKey imasskey(eInvMass, v_DetShapeType(), *this);
727 DNekMat &Imass = *GetStdMatrix(imasskey);
728
729 (*Mat) = Imass * Iprod;
730 }
731 break;
732 default:
733 {
735
736 if (mattype == eMass)
737 {
738 // For Fourier basis set the imaginary component
739 // of mean mode to have a unit diagonal component
740 // in mass matrix
742 {
743 (*Mat)(1, 1) = 1.0;
744 }
745 }
746 }
747 break;
748 }
749
750 return Mat;
751}
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

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

763{
764 if (outarray.size() != NumBndryCoeffs())
765 {
766 outarray = Array<OneD, unsigned int>(NumBndryCoeffs());
767 }
768 const LibUtilities::BasisType Btype = GetBasisType(0);
769 int nummodes = m_base[0]->GetNumModes();
770
771 outarray[0] = 0;
772
773 switch (Btype)
774 {
779 outarray[1] = nummodes - 1;
780 break;
783 outarray[1] = 1;
784 break;
785 default:
786 ASSERTL0(0, "Mapping array is not defined for this expansion");
787 break;
788 }
789}
@ 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.

Definition at line 630 of file StdSegExp.cpp.

633{
634 Blas::Dcopy(GetNumPoints(0), (m_base[0]->GetZ()).data(), 1, &coords_0[0],
635 1);
636}
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
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 907 of file StdSegExp.cpp.

912{
913 // parameters for higher dimnesion traces
914 if (maparray.size() != 1)
915 {
916 maparray = Array<OneD, unsigned int>(1);
917 }
918
919 maparray[0] = 0;
920
921 if (signarray.size() != 1)
922 {
923 signarray = Array<OneD, int>(1, 1);
924 }
925 else
926 {
927 signarray[0] = 1;
928 }
929}

◆ v_GetInteriorMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 791 of file StdSegExp.cpp.

792{
793 int i;
794 if (outarray.size() != GetNcoeffs() - NumBndryCoeffs())
795 {
796 outarray = Array<OneD, unsigned int>(GetNcoeffs() - NumBndryCoeffs());
797 }
798 const LibUtilities::BasisType Btype = GetBasisType(0);
799
800 switch (Btype)
801 {
806 for (i = 0; i < GetNcoeffs() - 2; i++)
807 {
808 outarray[i] = i + 1;
809 }
810 break;
813 for (i = 0; i < GetNcoeffs() - 2; i++)
814 {
815 outarray[i] = i + 2;
816 }
817 break;
818 default:
819 ASSERTL0(0, "Mapping array is not defined for this expansion");
820 break;
821 }
822}
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.

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

648{
649 return 2;
650}

◆ v_GetNverts()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 642 of file StdSegExp.cpp.

643{
644 return 2;
645}

◆ v_GetSimplexEquiSpacedConnectivity()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 841 of file StdSegExp.cpp.

843{
844 int np = m_base[0]->GetNumPoints();
845
846 conn = Array<OneD, int>(2 * (np - 1));
847 int cnt = 0;
848 for (int i = 0; i < np - 1; ++i)
849 {
850 conn[cnt++] = i;
851 conn[cnt++] = i + 1;
852 }
853}

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

861{
862 int order0 = m_base[0]->GetNumModes();
863
864 ASSERTL0(traceid < 2, "eid must be between 0 and 1");
865
866 if (maparray.size() != 1)
867 {
868 maparray = Array<OneD, unsigned int>(1);
869 }
870
871 const LibUtilities::BasisType bType = GetBasisType(0);
872
873 if (bType == LibUtilities::eModified_A)
874 {
875 maparray[0] = (traceid == 0) ? 0 : 1;
876 }
877 else if (bType == LibUtilities::eGLL_Lagrange ||
879 {
880 maparray[0] = (traceid == 0) ? 0 : order0 - 1;
881 }
882 else
883 {
884 ASSERTL0(false, "Unknown Basis");
885 }
886}

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

658{
659 return 0;
660}

◆ v_GetTraceNcoeffs()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 652 of file StdSegExp.cpp.

653{
654 return 1;
655}

◆ v_GetTraceNumPoints()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 662 of file StdSegExp.cpp.

663{
664 return 1;
665}

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

894{
895 v_GetTraceCoeffMap(tid, maparray);
896
897 if (signarray.size() != 1)
898 {
899 signarray = Array<OneD, int>(1, 1);
900 }
901 else
902 {
903 signarray[0] = 1;
904 }
905}
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.

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

826{
827 ASSERTL0((localVertexId == 0) || (localVertexId == 1),
828 "local vertex id"
829 "must be between 0 or 1");
830
831 int localDOF = localVertexId;
832
834 (localVertexId == 1))
835 {
836 localDOF = m_base[0]->GetNumModes() - 1;
837 }
838 return localDOF;
839}

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.

Definition at line 516 of file StdSegExp.cpp.

519{
520 int nquad = m_base[0]->GetNumPoints();
521
522 Array<OneD, NekDouble> physValues(nquad);
523 Array<OneD, NekDouble> dPhysValuesdx(nquad);
524 Array<OneD, NekDouble> wsp(m_ncoeffs);
525
526 v_BwdTrans(inarray, physValues);
527
528 // mass matrix operation
529 v_IProductWRTBase((m_base[0]->GetBdata()), physValues, wsp, 1);
530
531 // Laplacian matrix operation
532 v_PhysDeriv(physValues, dPhysValuesdx);
533 v_IProductWRTBase(m_base[0]->GetDbdata(), dPhysValuesdx, outarray, 1);
534 Blas::Daxpy(m_ncoeffs, mkey.GetConstFactor(eFactorLambda), wsp.data(), 1,
535 outarray.data(), 1);
536}
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.
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.

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::GetNumPoints(), Nektar::StdRegions::StdExpansion::m_base, Vmath::Vmul(), and Vmath::Vsum().

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

Definition at line 362 of file StdSegExp.cpp.

366{
367 int nquad = m_base[0]->GetNumPoints();
368 Array<OneD, NekDouble> tmp(nquad);
369 Array<OneD, const NekDouble> w = m_base[0]->GetW();
370
371 Vmath::Vmul(nquad, inarray, 1, w, 1, tmp, 1);
372
373 /* Comment below was a bug for collocated basis
374 if(coll_check&&m_base[0]->Collocation())
375 {
376 Vmath::Vcopy(nquad, tmp, 1, outarray, 1);
377 }
378 else
379 {
380 Blas::Dgemv('T',nquad,m_ncoeffs,1.0,base.data(),nquad,
381 &tmp[0],1,0.0,outarray.data(),1);
382 }*/
383
384 // Correct implementation
385 Blas::Dgemv('T', nquad, m_ncoeffs, 1.0, base.data(), nquad, &tmp[0], 1, 0.0,
386 outarray.data(), 1);
387}
std::vector< double > w(NPUPPER)

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

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

Definition at line 399 of file StdSegExp.cpp.

401{
402 v_IProductWRTBase(m_base[0]->GetBdata(), inarray, outarray, 1);
403}

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

423{
424 int nquad = m_base[0]->GetNumPoints();
425 Array<OneD, NekDouble> tmp(nquad);
426 Array<OneD, const NekDouble> w = m_base[0]->GetW();
427 Array<OneD, const NekDouble> base = m_base[0]->GetBdata();
428
429 if (multiplybyweights)
430 {
431 Vmath::Vmul(nquad, inarray, 1, w, 1, tmp, 1);
432
433 Blas::Dgemv('T', nquad, m_ncoeffs, 1.0, base.data(), nquad, &tmp[0], 1,
434 0.0, outarray.data(), 1);
435 }
436 else
437 {
438 Blas::Dgemv('T', nquad, m_ncoeffs, 1.0, base.data(), nquad, &inarray[0],
439 1, 0.0, outarray.data(), 1);
440 }
441}

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

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

Definition at line 405 of file StdSegExp.cpp.

408{
409 StdSegExp::IProductWRTDerivBase_SumFac(dir, inarray, outarray);
410}
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 412 of file StdSegExp.cpp.

415{
416 ASSERTL1(dir == 0, "input dir is out of range");
417 StdSegExp::v_IProductWRTBase(m_base[0]->GetDbdata(), inarray, outarray, 1);
418}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....

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() [1/2]

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.

Definition at line 492 of file StdSegExp.cpp.

495{
496 int nquad = m_base[0]->GetNumPoints();
497
498 Array<OneD, NekDouble> physValues(nquad);
499 Array<OneD, NekDouble> dPhysValuesdx(nquad);
500
501 v_BwdTrans(inarray, physValues);
502
503 // Laplacian matrix operation
504 v_PhysDeriv(physValues, dPhysValuesdx);
505 v_IProductWRTBase(m_base[0]->GetDbdata(), dPhysValuesdx, outarray, 1);
506}

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

◆ v_LaplacianMatrixOp() [2/2]

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 508 of file StdSegExp.cpp.

512{
513 StdExpansion::LaplacianMatrixOp_MatFree(k1, k2, inarray, outarray, mkey);
514}
void LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)

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

◆ v_LocCollapsedToLocCoord()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 452 of file StdSegExp.cpp.

454{
455 xi[0] = eta[0];
456}

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

448{
449 eta[0] = xi[0];
450}

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

622{
623 int nquad0 = m_base[0]->GetNumPoints();
624
625 const Array<OneD, const NekDouble> &w0 = m_base[0]->GetW();
626
627 Vmath::Vmul(nquad0, inarray.data(), 1, w0.data(), 1, outarray.data(), 1);
628}

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

◆ v_NumBndryCoeffs()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 667 of file StdSegExp.cpp.

668{
669 return 2;
670}

◆ v_NumDGBndryCoeffs()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 672 of file StdSegExp.cpp.

673{
674 return 2;
675}

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

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.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 475 of file StdSegExp.cpp.

479{
480 return StdExpansion1D::BaryTensorDeriv(coord, inarray, firstOrderDerivs);
481}
NekDouble BaryTensorDeriv(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)

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

◆ v_PhysEvalFirstSecondDeriv()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 483 of file StdSegExp.cpp.

488{
489 return StdExpansion1D::BaryTensorDeriv(coord, inarray, firstOrderDerivs,
490 secondOrderDerivs);
491}

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

471{
472 return StdExpansion::BaryEvaluateBasis<0>(coords[0], mode);
473}

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

189{
190 int n_coeffs = inarray.size();
191
192 Array<OneD, NekDouble> coeff(n_coeffs);
193 Array<OneD, NekDouble> coeff_tmp(n_coeffs, 0.0);
194 Array<OneD, NekDouble> tmp;
195 Array<OneD, NekDouble> tmp2;
196
197 int nmodes0 = m_base[0]->GetNumModes();
198
199 Vmath::Vcopy(n_coeffs, inarray, 1, coeff_tmp, 1);
200
201 const LibUtilities::PointsKey Pkey0(nmodes0,
203
204 LibUtilities::BasisKey b0(m_base[0]->GetBasisType(), nmodes0, Pkey0);
205
206 LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A, nmodes0, Pkey0);
207
208 LibUtilities::InterpCoeff1D(b0, coeff_tmp, bortho0, coeff);
209
210 Vmath::Zero(n_coeffs, coeff_tmp, 1);
211
212 Vmath::Vcopy(numMin, tmp = coeff, 1, tmp2 = coeff_tmp, 1);
213
214 LibUtilities::InterpCoeff1D(bortho0, coeff_tmp, b0, outarray);
215}
void InterpCoeff1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition PointsType.h:51
void Zero(int n, T *x, const int incx)
Zero vector.
Definition Vmath.hpp:273

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

◆ v_StdPhysDeriv()

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 538 of file StdSegExp.cpp.

540{
541 // Generate an orthogonal expansion
542 int nq = m_base[0]->GetNumPoints();
543 int nmodes = m_base[0]->GetNumModes();
544 // Declare orthogonal basis.
545 LibUtilities::PointsKey pKey(nq, m_base[0]->GetPointsType());
546
547 LibUtilities::BasisKey B(LibUtilities::eOrtho_A, nmodes, pKey);
548 StdSegExp OrthoExp(B);
549
550 // SVV parameters loaded from the .xml case file
551 NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDiffCoeff);
552 int cutoff = (int)(mkey.GetConstFactor(eFactorSVVCutoffRatio)) * nmodes;
553
554 Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
555
556 // project onto modal space.
557 OrthoExp.FwdTrans(array, orthocoeffs);
558
559 //
560 for (int j = 0; j < nmodes; ++j)
561 {
562 if (j >= cutoff) // to filter out only the "high-modes"
563 {
564 orthocoeffs[j] *=
565 (SvvDiffCoeff *
566 exp(-(j - nmodes) * (j - nmodes) /
567 ((NekDouble)((j - cutoff + 1) * (j - cutoff + 1)))));
568 }
569 else
570 {
571 orthocoeffs[j] *= 0.0;
572 }
573 }
574
575 // backward transform to physical space
576 OrthoExp.BwdTrans(orthocoeffs, array);
577}

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.