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

Class representing a segment element in reference space. More...

#include <StdSegExp.h>

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

Public Member Functions

 StdSegExp ()
 Default constructor. More...
 
 StdSegExp (const LibUtilities::BasisKey &Ba)
 Constructor using BasisKey class for quadrature points and order definition. More...
 
 StdSegExp (const StdSegExp &T)
 Copy Constructor. More...
 
 ~StdSegExp ()
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion1D
 StdExpansion1D ()
 
 StdExpansion1D (int numcoeffs, const LibUtilities::BasisKey &Ba)
 
 StdExpansion1D (const StdExpansion1D &T)
 
virtual ~StdExpansion1D ()
 
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...
 
- 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 (void) const
 
std::shared_ptr< StdExpansionGetLinStdExp (void) const
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion ()
 
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 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, 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 void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble >> &Fvec, Array< OneD, NekDouble > &outarray)
 
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual void v_DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
NekDouble Linf (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( L_\infty\) error \( |\epsilon|_\infty = \max |u - u_{exact}|\) where \( u_{exact}\) is given by the array sol. More...
 
NekDouble L2 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( L_2\) error, \( | \epsilon |_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 dx \right]^{1/2} d\xi_1 \) where \( u_{exact}\) is given by the array sol. More...
 
NekDouble H1 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( H^1\) error, \( | \epsilon |^1_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 + \nabla(u - u_{exact})\cdot\nabla(u - u_{exact})\cdot dx \right]^{1/2} d\xi_1 \) where \( u_{exact}\) is given by the array sol. More...
 
const LibUtilities::PointsKeyVector GetPointsKeys () const
 
DNekMatSharedPtr BuildInverseTransformationMatrix (const DNekScalMatSharedPtr &m_transformationmatrix)
 
void PhysInterpToSimplexEquiSpaced (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int npset=-1)
 This function performs an interpolation from the physical space points provided at input into an array of equispaced points which are not the collapsed coordinate. So for a tetrahedron you will only get a tetrahedral number of values. More...
 
void GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true)
 This function provides the connectivity of local simplices (triangles or tets) to connect the equispaced data points provided by PhysInterpToSimplexEquiSpaced. More...
 
void EquiSpacedToCoeffs (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs a projection/interpolation from the equispaced points sometimes used in post-processing onto the coefficient space. More...
 
template<class T >
std::shared_ptr< T > as ()
 
void IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
void GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 

Protected Member Functions

virtual NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray)
 Integrate the physical point list inarray over region and return the value. More...
 
virtual 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)
 Evaluate the derivative \( d/d{\xi_1} \) at the physical quadrature points given by inarray and return in outarray. More...
 
virtual void v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Calculate the derivative of the physical points in a given direction. More...
 
virtual 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)
 
virtual void v_StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Backward transform from coefficient space given in inarray and evaluate at the physical quadrature points outarray. More...
 
virtual void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Forward transform from physical quadrature space stored in inarray and evaluate the expansion coefficients and store in outarray. More...
 
virtual void v_BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_FwdTransBndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Inner product of inarray over region with respect to the expansion basis (this)->m_base[0] and return in outarray. More...
 
virtual void v_IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check)
 Inner product of inarray over region with respect to expansion basis base and return in outarray. More...
 
virtual void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
virtual void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 
virtual void v_LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi)
 
NekDouble v_PhysEvaluateBasis (const Array< OneD, const NekDouble > &coords, int mode) final
 
virtual void v_LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
virtual void v_ExponentialFilter (Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
 
virtual void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 
virtual void v_GetCoords (Array< OneD, NekDouble > &coords_0, Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2)
 
virtual void v_GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
virtual void v_GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
virtual int v_GetVertexMap (int localVertexId, bool useCoeffPacking=false)
 
virtual int v_GetNverts () const
 
virtual int v_GetNtraces () const
 
virtual int v_GetTraceNcoeffs (const int i) const
 
virtual int v_GetTraceNumPoints (const int i) const
 
virtual int v_NumBndryCoeffs () const
 
virtual int v_NumDGBndryCoeffs () const
 
virtual bool v_IsBoundaryInteriorExpansion ()
 
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
virtual LibUtilities::ShapeType v_DetShapeType () const
 Return Shape of region, using ShapeType enum list. i.e. Segment. More...
 
virtual DNekMatSharedPtr v_GenMatrix (const StdMatrixKey &mkey)
 
virtual DNekMatSharedPtr v_CreateStdMatrix (const StdMatrixKey &mkey)
 
virtual void v_GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true)
 
virtual void v_ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_GetTraceCoeffMap (const unsigned int traceid, Array< OneD, unsigned int > &maparray)
 Get the map of the coefficient location to teh local trace coefficients. More...
 
virtual void v_GetElmtTraceToTraceMap (const unsigned int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient, int P, int Q)
 
virtual void v_GetTraceToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient, int P, int Q)
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion1D
virtual NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals) override
 
- 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_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
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 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)
 
template<int DIR>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals)
 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)
 

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

Constructor & Destructor Documentation

◆ StdSegExp() [1/3]

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

Default constructor.

Definition at line 49 of file StdSegExp.cpp.

50 {
51 }

◆ StdSegExp() [2/3]

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

Constructor using BasisKey class for quadrature points and order definition.

Parameters
BaBasisKey class definition containing order and quadrature points.

Definition at line 60 of file StdSegExp.cpp.

61  : StdExpansion(Ba.GetNumModes(), 1, Ba),
62  StdExpansion1D(Ba.GetNumModes(), Ba)
63 {
64 }
StdExpansion()
Default Constructor.

◆ StdSegExp() [3/3]

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

Copy Constructor.

Definition at line 68 of file StdSegExp.cpp.

69 {
70 }

◆ ~StdSegExp()

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

Definition at line 72 of file StdSegExp.cpp.

73 {
74 }

Member Function Documentation

◆ v_BwdTrans()

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

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

207 {
208  int nquad = m_base[0]->GetNumPoints();
209 
210  if (m_base[0]->Collocation())
211  {
212  Vmath::Vcopy(nquad, inarray, 1, outarray, 1);
213  }
214  else
215  {
216 
217  Blas::Dgemv('N', nquad, m_base[0]->GetNumModes(), 1.0,
218  (m_base[0]->GetBdata()).get(), nquad, &inarray[0], 1, 0.0,
219  &outarray[0], 1);
220  }
221 }
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 = A x where A[m x n].
Definition: Blas.hpp:246
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 368 of file StdSegExp.cpp.

370 {
371  v_BwdTrans(inarray, outarray);
372 }
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Backward transform from coefficient space given in inarray and evaluate at the physical quadrature po...
Definition: StdSegExp.cpp:205

References v_BwdTrans().

◆ v_CalcNumberOfCoefficients()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 691 of file StdSegExp.cpp.

693 {
694  int nmodes = nummodes[modes_offset];
695  modes_offset += 1;
696 
697  return nmodes;
698 }

◆ v_CreateStdMatrix()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 767 of file StdSegExp.cpp.

768 {
769  return v_GenMatrix(mkey);
770 }
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey)
Definition: StdSegExp.cpp:704

References v_GenMatrix().

◆ v_DetShapeType()

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 79 of file StdSegExp.cpp.

80 {
82 }

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 596 of file StdSegExp.cpp.

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

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

◆ v_FillMode()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 498 of file StdSegExp.cpp.

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

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

◆ v_FwdTrans()

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

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

274 {
275  if (m_base[0]->Collocation())
276  {
277  Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
278  }
279  else
280  {
281  v_IProductWRTBase(inarray, outarray);
282 
283  // get Mass matrix inverse
284  StdMatrixKey masskey(eInvMass, v_DetShapeType(), *this);
285  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
286 
287  NekVector<NekDouble> in(m_ncoeffs, outarray, eCopy);
288  NekVector<NekDouble> out(m_ncoeffs, outarray, eWrapper);
289 
290  out = (*matsys) * in;
291  }
292 }
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:611
virtual LibUtilities::ShapeType v_DetShapeType() const
Return Shape of region, using ShapeType enum list. i.e. Segment.
Definition: StdSegExp.cpp:79
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Inner product of inarray over region with respect to the expansion basis (this)->m_base[0] and return...
Definition: StdSegExp.cpp:438
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 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 294 of file StdSegExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 704 of file StdSegExp.cpp.

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

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

Referenced by v_CreateStdMatrix().

◆ v_GetBoundaryMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 776 of file StdSegExp.cpp.

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

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

◆ v_GetCoords()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 647 of file StdSegExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 922 of file StdSegExp.cpp.

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

References Nektar::LibUtilities::P.

◆ v_GetInteriorMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 805 of file StdSegExp.cpp.

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

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

◆ v_GetNtraces()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 664 of file StdSegExp.cpp.

665 {
666  return 2;
667 }

◆ v_GetNverts()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 659 of file StdSegExp.cpp.

660 {
661  return 2;
662 }

◆ v_GetSimplexEquiSpacedConnectivity()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 855 of file StdSegExp.cpp.

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

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

◆ v_GetTraceCoeffMap()

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 874 of file StdSegExp.cpp.

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

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 669 of file StdSegExp.cpp.

670 {
671  boost::ignore_unused(i);
672  return 1;
673 }

◆ v_GetTraceNumPoints()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 675 of file StdSegExp.cpp.

676 {
677  boost::ignore_unused(i);
678  return 1;
679 }

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 903 of file StdSegExp.cpp.

907 {
908  boost::ignore_unused(orient, P, Q);
909 
910  v_GetTraceCoeffMap(tid, maparray);
911 
912  if (signarray.size() != 1)
913  {
914  signarray = Array<OneD, int>(1, 1);
915  }
916  else
917  {
918  signarray[0] = 1;
919  }
920 }
virtual void v_GetTraceCoeffMap(const unsigned int traceid, Array< OneD, unsigned int > &maparray)
Get the map of the coefficient location to teh local trace coefficients.
Definition: StdSegExp.cpp:874

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

◆ v_GetVertexMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 838 of file StdSegExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 533 of file StdSegExp.cpp.

536 {
537  int nquad = m_base[0]->GetNumPoints();
538 
539  Array<OneD, NekDouble> physValues(nquad);
540  Array<OneD, NekDouble> dPhysValuesdx(nquad);
541  Array<OneD, NekDouble> wsp(m_ncoeffs);
542 
543  v_BwdTrans(inarray, physValues);
544 
545  // mass matrix operation
546  v_IProductWRTBase((m_base[0]->GetBdata()), physValues, wsp, 1);
547 
548  // Laplacian matrix operation
549  v_PhysDeriv(physValues, dPhysValuesdx);
550  v_IProductWRTBase(m_base[0]->GetDbdata(), dPhysValuesdx, outarray, 1);
551  Blas::Daxpy(m_ncoeffs, mkey.GetConstFactor(eFactorLambda), wsp.get(), 1,
552  outarray.get(), 1);
553 }
virtual 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)
Evaluate the derivative at the physical quadrature points given by inarray and return in outarray.
Definition: StdSegExp.cpp:143
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:154

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

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

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

References 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 
)
protectedvirtual

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

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

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

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

440 {
441  v_IProductWRTBase(m_base[0]->GetBdata(), inarray, outarray, 1);
442 }

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

Referenced by v_FwdTrans(), v_FwdTransBndConstrained(), v_HelmholtzMatrixOp(), 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 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 460 of file StdSegExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 444 of file StdSegExp.cpp.

447 {
448  StdSegExp::IProductWRTDerivBase_SumFac(dir, inarray, outarray);
449 }
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 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 451 of file StdSegExp.cpp.

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

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

◆ v_IsBoundaryInteriorExpansion()

bool Nektar::StdRegions::StdSegExp::v_IsBoundaryInteriorExpansion ( )
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 84 of file StdSegExp.cpp.

85 {
86 
87  bool returnval = false;
88 
90  {
91  returnval = true;
92  }
93 
95  {
96  returnval = true;
97  }
98 
99  return returnval;
100 }

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 515 of file StdSegExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 492 of file StdSegExp.cpp.

494 {
495  xi[0] = eta[0];
496 }

◆ v_LocCoordToLocCollapsed()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 486 of file StdSegExp.cpp.

488 {
489  eta[0] = xi[0];
490 }

◆ v_MultiplyByStdQuadratureMetric()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 636 of file StdSegExp.cpp.

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

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

◆ v_NumBndryCoeffs()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 681 of file StdSegExp.cpp.

682 {
683  return 2;
684 }

◆ v_NumDGBndryCoeffs()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 686 of file StdSegExp.cpp.

687 {
688  return 2;
689 }

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

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

147 {
148  boost::ignore_unused(out_d1, out_d2);
149  PhysTensorDeriv(inarray, out_d0);
150 }
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 
)
protectedvirtual

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

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

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

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

511 {
512  return StdExpansion::BaryEvaluateBasis<0>(coords[0], mode);
513 }

◆ v_ReduceOrderCoeffs()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 223 of file StdSegExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 162 of file StdSegExp.cpp.

166 {
167  boost::ignore_unused(out_d1, out_d2);
168  PhysTensorDeriv(inarray, out_d0);
169  // PhysDeriv(inarray, out_d0);
170 }

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 172 of file StdSegExp.cpp.

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

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

◆ v_SVVLaplacianFilter()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 555 of file StdSegExp.cpp.

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

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.