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 FwdTrans_BndConstrained (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 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_FwdTrans_BndConstrained (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_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_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 50 of file StdSegExp.cpp.

51  {
52  }

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

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

◆ StdSegExp() [3/3]

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

Copy Constructor.

Definition at line 71 of file StdSegExp.cpp.

71  :
72  StdExpansion(T),
74  {
75  }

◆ ~StdSegExp()

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

Definition at line 78 of file StdSegExp.cpp.

79  {
80  }

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

226  {
227  int nquad = m_base[0]->GetNumPoints();
228 
229  if(m_base[0]->Collocation())
230  {
231  Vmath::Vcopy(nquad, inarray, 1, outarray, 1);
232  }
233  else
234  {
235 
236  Blas::Dgemv('N',nquad,m_base[0]->GetNumModes(),1.0,
237  (m_base[0]->GetBdata()).get(),
238  nquad,&inarray[0],1,0.0,&outarray[0],1);
239 
240  }
241  }
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:265
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199

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

395  {
396  v_BwdTrans(inarray, outarray);
397  }
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:223

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

735  {
736  int nmodes = nummodes[modes_offset];
737  modes_offset += 1;
738 
739  return nmodes;
740  }

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

813  {
814  return v_GenMatrix(mkey);
815  }
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey)
Definition: StdSegExp.cpp:746

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

86  {
88  }

References Nektar::LibUtilities::eSegment.

Referenced by v_FwdTrans(), v_FwdTrans_BndConstrained(), 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 629 of file StdSegExp.cpp.

634  {
635  // Generate an orthogonal expansion
636  int nq = m_base[0]->GetNumPoints();
637  int nmodes = m_base[0]->GetNumModes();
638  int P = nmodes - 1;
639  // Declare orthogonal basis.
640  LibUtilities::PointsKey pKey(nq,m_base[0]->GetPointsType());
641 
642  LibUtilities::BasisKey B(LibUtilities::eOrtho_A, nmodes, pKey);
643  StdSegExp OrthoExp(B);
644 
645  // Cutoff
646  int Pcut = cutoff*P;
647 
648  // Project onto orthogonal space.
649  Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
650  OrthoExp.FwdTrans(array,orthocoeffs);
651 
652  //
653  NekDouble fac;
654  for(int j = 0; j < nmodes; ++j)
655  {
656  //to filter out only the "high-modes"
657  if(j > Pcut)
658  {
659  fac = (NekDouble) (j - Pcut) / ( (NekDouble) (P - Pcut) );
660  fac = pow(fac, exponent);
661  orthocoeffs[j] *= exp(-alpha*fac);
662  }
663  }
664 
665  // backward transform to physical space
666  OrthoExp.BwdTrans(orthocoeffs,array);
667  }
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:208
StdSegExp()
Default constructor.
Definition: StdSegExp.cpp:50
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:45
double NekDouble
P
Definition: main.py:133

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 main::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 528 of file StdSegExp.cpp.

529  {
530  int nquad = m_base[0]->GetNumPoints();
531  const NekDouble * base = m_base[0]->GetBdata().get();
532 
533  ASSERTL2(mode <= m_ncoeffs,
534  "calling argument mode is larger than total expansion order");
535 
536  Vmath::Vcopy(nquad,(NekDouble *)base+mode*nquad,1, &outarray[0],1);
537  }
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:274

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

299  {
300  if(m_base[0]->Collocation())
301  {
302  Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
303  }
304  else
305  {
306  v_IProductWRTBase(inarray,outarray);
307 
308  // get Mass matrix inverse
309  StdMatrixKey masskey(eInvMass,v_DetShapeType(),*this);
310  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
311 
312  NekVector<NekDouble> in(m_ncoeffs,outarray,eCopy);
313  NekVector<NekDouble> out(m_ncoeffs,outarray,eWrapper);
314 
315  out = (*matsys)*in;
316  }
317  }
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:617
virtual LibUtilities::ShapeType v_DetShapeType() const
Return Shape of region, using ShapeType enum list. i.e. Segment.
Definition: StdSegExp.cpp:85
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:468
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69

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_FwdTrans_BndConstrained().

◆ v_FwdTrans_BndConstrained()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 319 of file StdSegExp.cpp.

322  {
323  if(m_base[0]->Collocation())
324  {
325  Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
326  }
327  else
328  {
329  int nInteriorDofs = m_ncoeffs-2;
330  int offset = 0;
331 
332  switch(m_base[0]->GetBasisType())
333  {
335  {
336  offset = 1;
337  }
338  break;
340  {
341  nInteriorDofs = m_ncoeffs;
342  offset = 0;
343  }
344  break;
347  {
348  offset = 2;
349  }
350  break;
351  default:
352  ASSERTL0(false,"This type of FwdTrans is not defined for this expansion type");
353  }
354 
355  fill(outarray.get(), outarray.get()+m_ncoeffs, 0.0 );
356 
358  {
359  outarray[GetVertexMap(0)] = inarray[0];
360  outarray[GetVertexMap(1)] = inarray[m_base[0]->GetNumPoints()-1];
361 
362  if(m_ncoeffs>2)
363  {
364  // ideally, we would like to have tmp0 to be replaced by
365  // outarray (currently MassMatrixOp does not allow aliasing)
366  Array<OneD, NekDouble> tmp0(m_ncoeffs);
367  Array<OneD, NekDouble> tmp1(m_ncoeffs);
368 
369  StdMatrixKey masskey(eMass,v_DetShapeType(),*this);
370  MassMatrixOp(outarray,tmp0,masskey);
371  v_IProductWRTBase(inarray,tmp1);
372 
373  Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
374 
375  // get Mass matrix inverse (only of interior DOF)
376  DNekMatSharedPtr matsys =
377  (m_stdStaticCondMatrixManager[masskey])-> GetBlock(1,1);
378 
379  Blas::Dgemv('N',nInteriorDofs,nInteriorDofs,1.0,
380  &(matsys->GetPtr())[0],nInteriorDofs,tmp1.get()
381  +offset,1,0.0,outarray.get()+offset,1);
382  }
383  }
384  else
385  {
386  StdSegExp::v_FwdTrans(inarray, outarray);
387  }
388  }
389 
390  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:158
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:762
int GetVertexMap(const int localVertexId, bool useCoeffPacking=false)
Definition: StdExpansion.h:697
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:297
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:49
@ eGauss_Lagrange
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:55
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:54
@ 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.cpp:372

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

747  {
748  DNekMatSharedPtr Mat;
749  MatrixType mattype;
750 
751  switch(mattype = mkey.GetMatrixType())
752  {
754  {
755  int nq = m_base[0]->GetNumPoints();
756 
757  // take definition from key
758  if(mkey.ConstFactorExists(eFactorConst))
759  {
760  nq = (int) mkey.GetConstFactor(eFactorConst);
761  }
762 
765  Array<OneD, NekDouble > coords (1);
766  DNekMatSharedPtr I ;
768  AllocateSharedPtr(neq, nq);
769 
770  for(int i = 0; i < neq; ++i)
771  {
772  coords[0] = -1.0 + 2*i/(NekDouble)(neq-1);
773  I = m_base[0]->GetI(coords);
774  Vmath::Vcopy(nq, I->GetRawPtr(), 1,
775  Mat->GetRawPtr()+i,neq);
776  }
777  }
778  break;
779  case eFwdTrans:
780  {
782  StdMatrixKey iprodkey(eIProductWRTBase,v_DetShapeType(),*this);
783  DNekMat &Iprod = *GetStdMatrix(iprodkey);
784  StdMatrixKey imasskey(eInvMass,v_DetShapeType(),*this);
785  DNekMat &Imass = *GetStdMatrix(imasskey);
786 
787  (*Mat) = Imass*Iprod;
788  }
789  break;
790  default:
791  {
793 
794  if(mattype == eMass)
795  {
796  // For Fourier basis set the imaginary component
797  // of mean mode to have a unit diagonal component
798  // in mass matrix
800  {
801  (*Mat)(1,1) = 1.0;
802  }
803  }
804  }
805  break;
806  }
807 
808  return Mat;
809  }
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:53
NekMatrix< NekDouble, StandardMatrixTag > DNekMat
Definition: NekTypeDefs.hpp:51

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

823  {
824  if(outarray.size() != NumBndryCoeffs())
825  {
826  outarray = Array<OneD, unsigned int>(NumBndryCoeffs());
827  }
828  const LibUtilities::BasisType Btype = GetBasisType(0);
829  int nummodes = m_base[0]->GetNumModes();
830 
831  outarray[0] = 0;
832 
833  switch(Btype)
834  {
839  outarray[1]= nummodes-1;
840  break;
843  outarray[1] = 1;
844  break;
845  default:
846  ASSERTL0(0,"Mapping array is not defined for this expansion");
847  break;
848  }
849  }
@ eChebyshev
Chebyshev Polynomials .
Definition: BasisType.h:57

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

686  {
687  boost::ignore_unused(coords_1, coords_2);
688  Blas::Dcopy(GetNumPoints(0),(m_base[0]->GetZ()).get(),
689  1,&coords_0[0],1);
690  }
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:221
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:160

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

◆ v_GetInteriorMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 851 of file StdSegExp.cpp.

852  {
853  int i;
854  if(outarray.size()!=GetNcoeffs()-NumBndryCoeffs())
855  {
856  outarray = Array<OneD, unsigned int>(GetNcoeffs()-NumBndryCoeffs());
857  }
858  const LibUtilities::BasisType Btype = GetBasisType(0);
859 
860  switch(Btype)
861  {
866  for(i = 0 ; i < GetNcoeffs()-2;i++)
867  {
868  outarray[i] = i+1;
869  }
870  break;
873  for(i = 0 ; i < GetNcoeffs()-2;i++)
874  {
875  outarray[i] = i+2;
876  }
877  break;
878  default:
879  ASSERTL0(0,"Mapping array is not defined for this expansion");
880  break;
881  }
882  }
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:124

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

◆ v_GetNtraces()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 705 of file StdSegExp.cpp.

706  {
707  return 2;
708  }

◆ v_GetNverts()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 700 of file StdSegExp.cpp.

701  {
702  return 2;
703  }

◆ v_GetSimplexEquiSpacedConnectivity()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 900 of file StdSegExp.cpp.

903  {
904  boost::ignore_unused(standard);
905  int np = m_base[0]->GetNumPoints();
906 
907  conn = Array<OneD, int>(2*(np-1));
908  int cnt = 0;
909  for(int i = 0; i < np-1; ++i)
910  {
911  conn[cnt++] = i;
912  conn[cnt++] = i+1;
913  }
914  }

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

◆ v_GetTraceNcoeffs()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 710 of file StdSegExp.cpp.

711  {
712  boost::ignore_unused(i);
713  return 1;
714  }

◆ v_GetTraceNumPoints()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 716 of file StdSegExp.cpp.

717  {
718  boost::ignore_unused(i);
719  return 1;
720  }

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

922  {
923  boost::ignore_unused(P,Q,orient);
924  int order0 = m_base[0]->GetNumModes();
925 
926  ASSERTL0(tid < 2,"eid must be between 0 and 1");
927 
928  if (maparray.size() != 2)
929  {
930  maparray = Array<OneD, unsigned int>(2);
931  }
932 
933  if(signarray.size() != 2)
934  {
935  signarray = Array<OneD, int>(2, 1);
936  }
937  else
938  {
939  fill(signarray.get(), signarray.get()+2, 1);
940  }
941 
942  const LibUtilities::BasisType bType = GetBasisType(0);
943 
944  if (bType == LibUtilities::eModified_A)
945  {
946  maparray[0] = 0;
947  maparray[1] = 1;
948  }
949  else if(bType == LibUtilities::eGLL_Lagrange ||
951  {
952  maparray[0] = 0;
953  maparray[1] = order0-1;
954  }
955  else
956  {
957  ASSERTL0(false,"Unknown Basis");
958  }
959  }

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

◆ v_GetVertexMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 884 of file StdSegExp.cpp.

885  {
886  boost::ignore_unused(useCoeffPacking);
887  ASSERTL0((localVertexId==0)||(localVertexId==1),"local vertex id"
888  "must be between 0 or 1");
889 
890  int localDOF = localVertexId;
891 
893  (localVertexId==1) )
894  {
895  localDOF = m_base[0]->GetNumModes()-1;
896  }
897  return localDOF;
898  }

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

570  {
571  int nquad = m_base[0]->GetNumPoints();
572 
573  Array<OneD,NekDouble> physValues(nquad);
574  Array<OneD,NekDouble> dPhysValuesdx(nquad);
575  Array<OneD,NekDouble> wsp(m_ncoeffs);
576 
577  v_BwdTrans(inarray,physValues);
578 
579  // mass matrix operation
580  v_IProductWRTBase((m_base[0]->GetBdata()),physValues,wsp,1);
581 
582  // Laplacian matrix operation
583  v_PhysDeriv(physValues,dPhysValuesdx);
584  v_IProductWRTBase(m_base[0]->GetDbdata(),dPhysValuesdx,outarray,1);
585  Blas::Daxpy(m_ncoeffs, mkey.GetConstFactor(eFactorLambda), wsp.get(), 1, outarray.get(), 1);
586  }
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:156
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:167

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

124  {
125  NekDouble Int = 0.0;
126  int nquad0 = m_base[0]->GetNumPoints();
127  Array<OneD, NekDouble> tmp(nquad0);
128  Array<OneD, const NekDouble> z = m_base[0]->GetZ();
129  Array<OneD, const NekDouble> w0 = m_base[0]->GetW();
130 
131  // multiply by integration constants
132  Vmath::Vmul(nquad0, inarray, 1, w0, 1, tmp, 1);
133 
134  Int = Vmath::Vsum(nquad0, tmp, 1);
135 
136  return Int;
137  }
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:192
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:846

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

433  {
434  boost::ignore_unused(coll_check);
435 
436  int nquad = m_base[0]->GetNumPoints();
437  Array<OneD, NekDouble> tmp(nquad);
438  Array<OneD, const NekDouble> w = m_base[0]->GetW();
439 
440  Vmath::Vmul(nquad, inarray, 1, w, 1, tmp, 1);
441 
442  /* Comment below was a bug for collocated basis
443  if(coll_check&&m_base[0]->Collocation())
444  {
445  Vmath::Vcopy(nquad, tmp, 1, outarray, 1);
446  }
447  else
448  {
449  Blas::Dgemv('T',nquad,m_ncoeffs,1.0,base.get(),nquad,
450  &tmp[0],1,0.0,outarray.get(),1);
451  }*/
452 
453  // Correct implementation
454  Blas::Dgemv('T',nquad,m_ncoeffs,1.0,base.get(),nquad,
455  &tmp[0],1,0.0,outarray.get(),1);
456  }

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

470  {
471  v_IProductWRTBase(m_base[0]->GetBdata(),inarray,outarray,1);
472  }

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

Referenced by v_FwdTrans(), v_FwdTrans_BndConstrained(), v_HelmholtzMatrixOp(), v_IProductWRTDerivBase(), 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 484 of file StdSegExp.cpp.

488  {
489  int nquad = m_base[0]->GetNumPoints();
490  Array<OneD, NekDouble> tmp(nquad);
491  Array<OneD, const NekDouble> w = m_base[0]->GetW();
492  Array<OneD, const NekDouble> base = m_base[0]->GetBdata();
493 
494  if(multiplybyweights)
495  {
496  Vmath::Vmul(nquad, inarray, 1, w, 1, tmp, 1);
497 
498  Blas::Dgemv('T',nquad,m_ncoeffs,1.0,base.get(),nquad,
499  &tmp[0],1,0.0,outarray.get(),1);
500  }
501  else
502  {
503  Blas::Dgemv('T',nquad,m_ncoeffs,1.0,base.get(),nquad,
504  &inarray[0],1,0.0,outarray.get(),1);
505  }
506  }

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

478  {
479  boost::ignore_unused(dir);
480  ASSERTL1(dir >= 0 && dir < 1,"input dir is out of range");
481  v_IProductWRTBase(m_base[0]->GetDbdata(),inarray,outarray,1);
482  }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250

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

91  {
92 
93  bool returnval = false;
94 
96  {
97  returnval = true;
98  }
99 
101  {
102  returnval = true;
103  }
104 
105  return returnval;
106  }

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

550  {
551  boost::ignore_unused(mkey);
552 
553  int nquad = m_base[0]->GetNumPoints();
554 
555  Array<OneD,NekDouble> physValues(nquad);
556  Array<OneD,NekDouble> dPhysValuesdx(nquad);
557 
558  v_BwdTrans(inarray,physValues);
559 
560  // Laplacian matrix operation
561  v_PhysDeriv(physValues,dPhysValuesdx);
562  v_IProductWRTBase(m_base[0]->GetDbdata(),dPhysValuesdx,outarray,1);
563  }

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

524  {
525  xi[0] = eta[0];
526  }

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

517  {
518  eta[0] = xi[0];
519  }

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

673  {
674  int nquad0 = m_base[0]->GetNumPoints();
675 
676  const Array<OneD, const NekDouble>& w0 = m_base[0]->GetW();
677 
678  Vmath::Vmul(nquad0, inarray.get(),1,
679  w0.get(),1,outarray.get(),1);
680  }

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

723  {
724  return 2;
725  }

◆ v_NumDGBndryCoeffs()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 727 of file StdSegExp.cpp.

728  {
729  return 2;
730  }

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

160  {
161  boost::ignore_unused(out_d1, out_d2);
162  PhysTensorDeriv(inarray,out_d0);
163  }
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 166 of file StdSegExp.cpp.

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

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

542  {
543  return StdExpansion::BaryEvaluateBasis<0>(coords[0], mode);
544  }

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

247  {
248  int n_coeffs = inarray.size();
249 
250  Array<OneD, NekDouble> coeff(n_coeffs);
251  Array<OneD, NekDouble> coeff_tmp(n_coeffs,0.0);
252  Array<OneD, NekDouble> tmp;
253  Array<OneD, NekDouble> tmp2;
254 
255  int nmodes0 = m_base[0]->GetNumModes();
256 
257  Vmath::Vcopy(n_coeffs,inarray,1,coeff_tmp,1);
258 
259  const LibUtilities::PointsKey Pkey0(
261 
262  LibUtilities::BasisKey b0(m_base[0]->GetBasisType(),nmodes0,Pkey0);
263 
264  LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A,nmodes0,Pkey0);
265 
267  b0, coeff_tmp, bortho0, coeff);
268 
269  Vmath::Zero(n_coeffs,coeff_tmp,1);
270 
271  Vmath::Vcopy(numMin,
272  tmp = coeff,1,
273  tmp2 = coeff_tmp,1);
274 
276  bortho0, coeff_tmp, b0, outarray);
277  }
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:51
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:436

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

181  {
182  boost::ignore_unused(out_d1, out_d2);
183  PhysTensorDeriv(inarray,out_d0);
184  // PhysDeriv(inarray, out_d0);
185  }

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

191  {
192  boost::ignore_unused(dir);
193  ASSERTL1(dir==0,"input dir is out of range");
194  PhysTensorDeriv(inarray,outarray);
195  // PhysDeriv(inarray, outarray);
196  }

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

590  {
591  // Generate an orthogonal expansion
592  int nq = m_base[0]->GetNumPoints();
593  int nmodes = m_base[0]->GetNumModes();
594  // Declare orthogonal basis.
595  LibUtilities::PointsKey pKey(nq,m_base[0]->GetPointsType());
596 
597  LibUtilities::BasisKey B(LibUtilities::eOrtho_A, nmodes, pKey);
598  StdSegExp OrthoExp(B);
599 
600  //SVV parameters loaded from the .xml case file
601  NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDiffCoeff);
602  int cutoff = (int) (mkey.GetConstFactor(eFactorSVVCutoffRatio))*nmodes;
603 
604  Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
605 
606  // project onto modal space.
607  OrthoExp.FwdTrans(array,orthocoeffs);
608 
609  //
610  for(int j = 0; j < nmodes; ++j)
611  {
612  if(j >= cutoff)//to filter out only the "high-modes"
613  {
614  orthocoeffs[j] *=
615  (SvvDiffCoeff*exp(-(j-nmodes)*(j-nmodes)/
616  ((NekDouble)((j-cutoff+1)*
617  (j-cutoff+1)))));
618  }
619  else
620  {
621  orthocoeffs[j] *= 0.0;
622  }
623  }
624 
625  // backward transform to physical space
626  OrthoExp.BwdTrans(orthocoeffs,array);
627  }

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.