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

#include <StdTriExp.h>

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

Public Member Functions

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

Protected Member Functions

virtual NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray) override
 Integrates the specified function over the domain. More...
 
void v_PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray) override
 Calculate the derivative of the physical points. More...
 
virtual void v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Calculate the derivative of the physical points in a given direction. More...
 
virtual void v_StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray) override
 
virtual void v_StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Backward tranform for triangular elements. More...
 
virtual void v_BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
virtual void v_BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1) override
 
virtual void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Transform a given function from physical quadrature space to coefficient space. More...
 
virtual void v_FwdTransBndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Calculate the inner product of inarray with respect to the basis B=base0[p]*base1[pq] and put into outarray. More...
 
virtual void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
 
virtual void v_IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1) override
 
virtual void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
virtual void v_IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
virtual void v_LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta) override
 
virtual void v_LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi) override
 
virtual void v_FillMode (const int mode, Array< OneD, NekDouble > &outarray) override
 
NekDouble v_PhysEvaluateBasis (const Array< OneD, const NekDouble > &coords, int mode) final override
 
NekDouble v_PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
 
virtual int v_GetNverts () const final override
 
virtual int v_GetNtraces () const final override
 
virtual LibUtilities::ShapeType v_DetShapeType () const final override
 
virtual int v_NumBndryCoeffs () const override
 
virtual int v_NumDGBndryCoeffs () const override
 
virtual int v_GetTraceNcoeffs (const int i) const override
 
virtual int v_GetTraceIntNcoeffs (const int i) const override
 
virtual int v_GetTraceNumPoints (const int i) const override
 
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset) override
 
virtual void v_GetCoords (Array< OneD, NekDouble > &coords_x, Array< OneD, NekDouble > &coords_y, Array< OneD, NekDouble > &coords_z) override
 
virtual bool v_IsBoundaryInteriorExpansion () const override
 
virtual const LibUtilities::BasisKey v_GetTraceBasisKey (const int i, const int j) const override
 
virtual int v_GetVertexMap (int localVertexId, bool useCoeffPacking=false) override
 
virtual void v_GetInteriorMap (Array< OneD, unsigned int > &outarray) override
 
virtual void v_GetBoundaryMap (Array< OneD, unsigned int > &outarray) override
 
virtual void v_GetTraceCoeffMap (const unsigned int traceid, Array< OneD, unsigned int > &maparray) override
 
void v_GetTraceInteriorToElementMap (const int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation edgeOrient=eForwards) override
 
virtual DNekMatSharedPtr v_GenMatrix (const StdMatrixKey &mkey) override
 
virtual DNekMatSharedPtr v_CreateStdMatrix (const StdMatrixKey &mkey) override
 
void v_MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
 
virtual void v_LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
 
virtual void v_SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey) override
 
virtual void v_ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
virtual void v_LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
 
virtual void v_WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
 
virtual void v_HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
 
virtual void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
virtual void v_GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true) override
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion2D
virtual NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals) override
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
virtual NekDouble v_PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals) override
 
virtual void v_LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
virtual void v_HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
virtual void v_GetElmtTraceToTraceMap (const unsigned int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient, int P, int Q) override
 Determine the mapping to re-orientate the coefficients along the element trace (assumed to align with the standard element) into the orientation of the local trace given by edgeOrient. More...
 
virtual void v_GetTraceToElementMap (const int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient=eForwards, int P=-1, int Q=-1) override
 
virtual void v_GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat) override
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion
DNekMatSharedPtr CreateStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr CreateStdStaticCondMatrix (const StdMatrixKey &mkey)
 Create the static condensation of a matrix when using a boundary interior decomposition. More...
 
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase_SumFac (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GeneralMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
 
void LaplacianMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp_MatFree (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void 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)
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals, NekDouble &deriv, NekDouble &deriv2)
 This function performs the barycentric interpolation of the polynomial stored in coord at a point physvals using barycentric interpolation weights in direction. More...
 
template<int DIR>
NekDouble BaryEvaluateBasis (const NekDouble &coord, const int &mode)
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals)
 Helper function to pass an unused value by reference into BaryEvaluate. More...
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals, NekDouble &deriv)
 

Additional Inherited Members

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

Detailed Description

Definition at line 48 of file StdTriExp.h.

Constructor & Destructor Documentation

◆ StdTriExp() [1/3]

Nektar::StdRegions::StdTriExp::StdTriExp ( )

Definition at line 48 of file StdTriExp.cpp.

49 {
50 }

◆ StdTriExp() [2/3]

Nektar::StdRegions::StdTriExp::StdTriExp ( const LibUtilities::BasisKey Ba,
const LibUtilities::BasisKey Bb 
)

Definition at line 52 of file StdTriExp.cpp.

55  Ba.GetNumModes(), Bb.GetNumModes()),
56  2, Ba, Bb),
58  Ba.GetNumModes(), Bb.GetNumModes()),
59  Ba, Bb)
60 {
61  ASSERTL0(Ba.GetNumModes() <= Bb.GetNumModes(),
62  "order in 'a' direction is higher than order "
63  "in 'b' direction");
64 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
StdExpansion()
Default Constructor.
int getNumberOfCoefficients(int Na, int Nb)
Definition: ShapeType.hpp:114

References ASSERTL0, and Nektar::LibUtilities::BasisKey::GetNumModes().

◆ StdTriExp() [3/3]

Nektar::StdRegions::StdTriExp::StdTriExp ( const StdTriExp T)

Definition at line 66 of file StdTriExp.cpp.

67 {
68 }

◆ ~StdTriExp()

Nektar::StdRegions::StdTriExp::~StdTriExp ( )
overridevirtual

Definition at line 71 of file StdTriExp.cpp.

72 {
73 }

Member Function Documentation

◆ v_BwdTrans()

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

Backward tranform for triangular elements.

Note
'q' (base[1]) runs fastest in this element.

Implements Nektar::StdRegions::StdExpansion.

Definition at line 240 of file StdTriExp.cpp.

242 {
243  v_BwdTrans_SumFac(inarray, outarray);
244 }
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdTriExp.cpp:246

References v_BwdTrans_SumFac().

◆ v_BwdTrans_SumFac()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp, and Nektar::LocalRegions::NodalTriExp.

Definition at line 246 of file StdTriExp.cpp.

248 {
249  Array<OneD, NekDouble> wsp(m_base[0]->GetNumPoints() *
250  m_base[1]->GetNumModes());
251 
252  BwdTrans_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetBdata(), inarray,
253  outarray, wsp);
254 }
void BwdTrans_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:224
Array< OneD, LibUtilities::BasisSharedPtr > m_base

References Nektar::StdRegions::StdExpansion2D::BwdTrans_SumFacKernel(), Nektar::StdRegions::StdExpansion::GetNumPoints(), and Nektar::StdRegions::StdExpansion::m_base.

Referenced by v_BwdTrans(), and Nektar::StdRegions::StdNodalTriExp::v_BwdTrans_SumFac().

◆ v_BwdTrans_SumFacKernel()

void Nektar::StdRegions::StdTriExp::v_BwdTrans_SumFacKernel ( const Array< OneD, const NekDouble > &  base0,
const Array< OneD, const NekDouble > &  base1,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
Array< OneD, NekDouble > &  wsp,
bool  doCheckCollDir0,
bool  doCheckCollDir1 
)
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion2D.

Definition at line 256 of file StdTriExp.cpp.

262 {
263  boost::ignore_unused(doCheckCollDir0, doCheckCollDir1);
264 
265  int i;
266  int mode;
267  int nquad0 = m_base[0]->GetNumPoints();
268  int nquad1 = m_base[1]->GetNumPoints();
269  int nmodes0 = m_base[0]->GetNumModes();
270  int nmodes1 = m_base[1]->GetNumModes();
271 
272  ASSERTL1(wsp.size() >= nquad0 * nmodes1,
273  "Workspace size is not sufficient");
276  "Basis[1] is not of general tensor type");
277 
278  for (i = mode = 0; i < nmodes0; ++i)
279  {
280  Blas::Dgemv('N', nquad1, nmodes1 - i, 1.0, base1.get() + mode * nquad1,
281  nquad1, &inarray[0] + mode, 1, 0.0, &wsp[0] + i * nquad1,
282  1);
283  mode += nmodes1 - i;
284  }
285 
286  // fix for modified basis by splitting top vertex mode
288  {
289  Blas::Daxpy(nquad1, inarray[1], base1.get() + nquad1, 1,
290  &wsp[0] + nquad1, 1);
291  }
292 
293  Blas::Dgemm('N', 'T', nquad0, nquad1, nmodes0, 1.0, base0.get(), nquad0,
294  &wsp[0], nquad1, 0.0, &outarray[0], nquad0);
295 }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:272
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:162
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
static void Dgemm(const char &transa, const char &transb, const int &m, const int &n, const int &k, const double &alpha, const double *a, const int &lda, const double *b, const int &ldb, const double &beta, double *c, const int &ldc)
BLAS level 3: Matrix-matrix multiply C = A x B where op(A)[m x k], op(B)[k x n], C[m x n] DGEMM perfo...
Definition: Blas.hpp:368
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
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:51
@ eOrtho_B
Principle Orthogonal Functions .
Definition: BasisType.h:46
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:50

References ASSERTL1, ASSERTL2, Blas::Daxpy(), Blas::Dgemm(), Blas::Dgemv(), Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::LibUtilities::eOrtho_B, Nektar::StdRegions::StdExpansion::GetBasisType(), and Nektar::StdRegions::StdExpansion::m_base.

◆ v_CalcNumberOfCoefficients()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 819 of file StdTriExp.cpp.

821 {
823  nummodes[modes_offset], nummodes[modes_offset + 1]);
824  modes_offset += 2;
825 
826  return nmodes;
827 }

References Nektar::LibUtilities::StdTriData::getNumberOfCoefficients().

◆ v_CreateStdMatrix()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp, Nektar::LocalRegions::NodalTriExp, and Nektar::StdRegions::StdNodalTriExp.

Definition at line 1281 of file StdTriExp.cpp.

1282 {
1283  return v_GenMatrix(mkey);
1284 }
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override
Definition: StdTriExp.cpp:1207

References v_GenMatrix().

◆ v_DetShapeType()

LibUtilities::ShapeType Nektar::StdRegions::StdTriExp::v_DetShapeType ( ) const
finaloverrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 752 of file StdTriExp.cpp.

753 {
755 }

References Nektar::LibUtilities::eTriangle.

◆ v_FillMode()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp.

Definition at line 628 of file StdTriExp.cpp.

629 {
630  int i, m;
631  int nquad0 = m_base[0]->GetNumPoints();
632  int nquad1 = m_base[1]->GetNumPoints();
633  int order0 = m_base[0]->GetNumModes();
634  int order1 = m_base[1]->GetNumModes();
635  int mode0 = 0;
636  Array<OneD, const NekDouble> base0 = m_base[0]->GetBdata();
637  Array<OneD, const NekDouble> base1 = m_base[1]->GetBdata();
638 
639  ASSERTL2(mode <= m_ncoeffs, "calling argument mode is larger than "
640  "total expansion order");
641 
642  m = order1;
643  for (i = 0; i < order0; ++i, m += order1 - i)
644  {
645  if (m > mode)
646  {
647  mode0 = i;
648  break;
649  }
650  }
651 
652  // deal with top vertex mode in modified basis
653  if (mode == 1 && m_base[0]->GetBasisType() == LibUtilities::eModified_A)
654  {
655  Vmath::Fill(nquad0 * nquad1, 1.0, outarray, 1);
656  }
657  else
658  {
659  for (i = 0; i < nquad1; ++i)
660  {
661  Vmath::Vcopy(nquad0, (NekDouble *)(base0.get() + mode0 * nquad0), 1,
662  &outarray[0] + i * nquad0, 1);
663  }
664  }
665 
666  for (i = 0; i < nquad0; ++i)
667  {
668  Vmath::Vmul(nquad1, (NekDouble *)(base1.get() + mode * nquad1), 1,
669  &outarray[0] + i, nquad0, &outarray[0] + i, nquad0);
670  }
671 }
double NekDouble
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
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255

References ASSERTL2, Nektar::LibUtilities::eModified_A, Vmath::Fill(), Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Vmath::Vcopy(), and Vmath::Vmul().

Referenced by Nektar::StdRegions::StdNodalTriExp::GenNBasisTransMatrix().

◆ v_FwdTrans()

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

Transform a given function from physical quadrature space to coefficient space.

See also
StdExpansion::FwdTrans

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp, Nektar::LocalRegions::TriExp, and Nektar::LocalRegions::NodalTriExp.

Definition at line 297 of file StdTriExp.cpp.

299 {
300  v_IProductWRTBase(inarray, outarray);
301 
302  // get Mass matrix inverse
303  StdMatrixKey masskey(eInvMass, DetShapeType(), *this);
304  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
305 
306  // copy inarray in case inarray == outarray
307  NekVector<NekDouble> in(m_ncoeffs, outarray, eCopy);
308  NekVector<NekDouble> out(m_ncoeffs, outarray, eWrapper);
309 
310  out = (*matsys) * in;
311 }
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:609
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:373
void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Calculate the inner product of inarray with respect to the basis B=base0[p]*base1[pq] and put into ou...
Definition: StdTriExp.cpp:444
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75

References Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::eCopy, Nektar::StdRegions::eInvMass, Nektar::eWrapper, Nektar::StdRegions::StdExpansion::GetStdMatrix(), Nektar::StdRegions::StdExpansion::m_ncoeffs, and v_IProductWRTBase().

◆ v_FwdTransBndConstrained()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp.

Definition at line 313 of file StdTriExp.cpp.

316 {
317  int i, j;
318  int npoints[2] = {m_base[0]->GetNumPoints(), m_base[1]->GetNumPoints()};
319  int nmodes[2] = {m_base[0]->GetNumModes(), m_base[1]->GetNumModes()};
320 
321  fill(outarray.get(), outarray.get() + m_ncoeffs, 0.0);
322 
323  Array<OneD, NekDouble> physEdge[3];
324  Array<OneD, NekDouble> coeffEdge[3];
325  for (i = 0; i < 3; i++)
326  {
327  physEdge[i] = Array<OneD, NekDouble>(npoints[i != 0]);
328  coeffEdge[i] = Array<OneD, NekDouble>(nmodes[i != 0]);
329  }
330 
331  for (i = 0; i < npoints[0]; i++)
332  {
333  physEdge[0][i] = inarray[i];
334  }
335 
336  for (i = 0; i < npoints[1]; i++)
337  {
338  physEdge[1][i] = inarray[npoints[0] - 1 + i * npoints[0]];
339  physEdge[2][i] =
340  inarray[(npoints[1] - 1) * npoints[0] - i * npoints[0]];
341  }
342 
343  StdSegExpSharedPtr segexp[2] = {
345  m_base[0]->GetBasisKey()),
347  m_base[1]->GetBasisKey())};
348 
349  Array<OneD, unsigned int> mapArray;
350  Array<OneD, int> signArray;
351  NekDouble sign;
352 
353  for (i = 0; i < 3; i++)
354  {
355  segexp[i != 0]->FwdTransBndConstrained(physEdge[i], coeffEdge[i]);
356 
357  GetTraceToElementMap(i, mapArray, signArray);
358  for (j = 0; j < nmodes[i != 0]; j++)
359  {
360  sign = (NekDouble)signArray[j];
361  outarray[mapArray[j]] = sign * coeffEdge[i][j];
362  }
363  }
364 
365  Array<OneD, NekDouble> tmp0(m_ncoeffs);
366  Array<OneD, NekDouble> tmp1(m_ncoeffs);
367 
368  StdMatrixKey masskey(eMass, DetShapeType(), *this);
369  MassMatrixOp(outarray, tmp0, masskey);
370  v_IProductWRTBase(inarray, tmp1);
371 
372  Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
373 
374  // get Mass matrix inverse (only of interior DOF)
375  // use block (1,1) of the static condensed system
376  // note: this block alreay contains the inverse matrix
377  DNekMatSharedPtr matsys =
378  (m_stdStaticCondMatrixManager[masskey])->GetBlock(1, 1);
379 
380  int nBoundaryDofs = v_NumBndryCoeffs();
381  int nInteriorDofs = m_ncoeffs - nBoundaryDofs;
382 
383  Array<OneD, NekDouble> rhs(nInteriorDofs);
384  Array<OneD, NekDouble> result(nInteriorDofs);
385 
386  v_GetInteriorMap(mapArray);
387 
388  for (i = 0; i < nInteriorDofs; i++)
389  {
390  rhs[i] = tmp1[mapArray[i]];
391  }
392 
393  Blas::Dgemv('N', nInteriorDofs, nInteriorDofs, 1.0, &(matsys->GetPtr())[0],
394  nInteriorDofs, rhs.get(), 1, 0.0, result.get(), 1);
395 
396  for (i = 0; i < nInteriorDofs; i++)
397  {
398  outarray[mapArray[i]] = result[i];
399  }
400 }
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:49
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:758
void GetTraceToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
Definition: StdExpansion.h:690
LibUtilities::NekManager< StdMatrixKey, DNekBlkMat, StdMatrixKey::opLess > m_stdStaticCondMatrixManager
virtual void v_GetInteriorMap(Array< OneD, unsigned int > &outarray) override
Definition: StdTriExp.cpp:1018
virtual int v_NumBndryCoeffs() const override
Definition: StdTriExp.cpp:757
std::shared_ptr< StdSegExp > StdSegExpSharedPtr
Definition: StdSegExp.h:267
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 Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::StdRegions::StdExpansion::DetShapeType(), Blas::Dgemv(), Nektar::StdRegions::eMass, Nektar::StdRegions::StdExpansion::GetTraceToElementMap(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::StdRegions::StdExpansion::m_stdStaticCondMatrixManager, Nektar::StdRegions::StdExpansion::MassMatrixOp(), sign, v_GetInteriorMap(), v_IProductWRTBase(), v_NumBndryCoeffs(), and Vmath::Vsub().

◆ v_GenMatrix()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp, Nektar::LocalRegions::NodalTriExp, and Nektar::StdRegions::StdNodalTriExp.

Definition at line 1207 of file StdTriExp.cpp.

1208 {
1209 
1210  MatrixType mtype = mkey.GetMatrixType();
1211 
1212  DNekMatSharedPtr Mat;
1213 
1214  switch (mtype)
1215  {
1217  {
1218  int nq0, nq1, nq;
1219 
1220  nq0 = m_base[0]->GetNumPoints();
1221  nq1 = m_base[1]->GetNumPoints();
1222 
1223  // take definition from key
1224  if (mkey.ConstFactorExists(eFactorConst))
1225  {
1226  nq = (int)mkey.GetConstFactor(eFactorConst);
1227  }
1228  else
1229  {
1230  nq = max(nq0, nq1);
1231  }
1232 
1234  Array<OneD, Array<OneD, NekDouble>> coords(neq);
1235  Array<OneD, NekDouble> coll(2);
1236  Array<OneD, DNekMatSharedPtr> I(2);
1237  Array<OneD, NekDouble> tmp(nq0);
1238 
1239  Mat = MemoryManager<DNekMat>::AllocateSharedPtr(neq, nq0 * nq1);
1240  int cnt = 0;
1241 
1242  for (int i = 0; i < nq; ++i)
1243  {
1244  for (int j = 0; j < nq - i; ++j, ++cnt)
1245  {
1246  coords[cnt] = Array<OneD, NekDouble>(2);
1247  coords[cnt][0] = -1.0 + 2 * j / (NekDouble)(nq - 1);
1248  coords[cnt][1] = -1.0 + 2 * i / (NekDouble)(nq - 1);
1249  }
1250  }
1251 
1252  for (int i = 0; i < neq; ++i)
1253  {
1254  LocCoordToLocCollapsed(coords[i], coll);
1255 
1256  I[0] = m_base[0]->GetI(coll);
1257  I[1] = m_base[1]->GetI(coll + 1);
1258 
1259  // interpolate first coordinate direction
1260  for (int j = 0; j < nq1; ++j)
1261  {
1262  NekDouble fac = (I[1]->GetPtr())[j];
1263  Vmath::Smul(nq0, fac, I[0]->GetPtr(), 1, tmp, 1);
1264 
1265  Vmath::Vcopy(nq0, &tmp[0], 1,
1266  Mat->GetRawPtr() + j * nq0 * neq + i, neq);
1267  }
1268  }
1269  break;
1270  }
1271  default:
1272  {
1274  break;
1275  }
1276  }
1277 
1278  return Mat;
1279 }
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:248

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::StdRegions::StdMatrixKey::ConstFactorExists(), Nektar::StdRegions::StdExpansion::CreateGeneralMatrix(), Nektar::StdRegions::eFactorConst, Nektar::StdRegions::ePhysInterpToEquiSpaced, Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdMatrixKey::GetMatrixType(), Nektar::LibUtilities::StdTriData::getNumberOfCoefficients(), Nektar::StdRegions::StdExpansion::LocCoordToLocCollapsed(), Nektar::StdRegions::StdExpansion::m_base, Vmath::Smul(), and Vmath::Vcopy().

Referenced by v_CreateStdMatrix().

◆ v_GetBoundaryMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp.

Definition at line 1048 of file StdTriExp.cpp.

1049 {
1052  "Expansion not of expected type");
1053  int i;
1054  int cnt;
1055  int nummodes0, nummodes1;
1056  int value;
1057 
1058  if (outarray.size() != NumBndryCoeffs())
1059  {
1060  outarray = Array<OneD, unsigned int>(NumBndryCoeffs());
1061  }
1062 
1063  nummodes0 = m_base[0]->GetNumModes();
1064  nummodes1 = m_base[1]->GetNumModes();
1065 
1066  value = 2 * nummodes1 - 1;
1067  for (i = 0; i < value; i++)
1068  {
1069  outarray[i] = i;
1070  }
1071  cnt = value;
1072 
1073  for (i = 0; i < nummodes0 - 2; i++)
1074  {
1075  outarray[cnt++] = value;
1076  value += nummodes1 - 2 - i;
1077  }
1078 }

References ASSERTL1, 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::StdTriExp::v_GetCoords ( Array< OneD, NekDouble > &  coords_x,
Array< OneD, NekDouble > &  coords_y,
Array< OneD, NekDouble > &  coords_z 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp, and Nektar::LocalRegions::NodalTriExp.

Definition at line 829 of file StdTriExp.cpp.

832 {
833  boost::ignore_unused(coords_2);
834 
835  Array<OneD, const NekDouble> z0 = m_base[0]->GetZ();
836  Array<OneD, const NekDouble> z1 = m_base[1]->GetZ();
837  int nq0 = GetNumPoints(0);
838  int nq1 = GetNumPoints(1);
839  int i, j;
840 
841  for (i = 0; i < nq1; ++i)
842  {
843  for (j = 0; j < nq0; ++j)
844  {
845  coords_0[i * nq0 + j] = (1 + z0[j]) * (1 - z1[i]) / 2.0 - 1.0;
846  }
847  Vmath::Fill(nq0, z1[i], &coords_1[0] + i * nq0, 1);
848  }
849 }

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

◆ v_GetInteriorMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp.

Definition at line 1018 of file StdTriExp.cpp.

1019 {
1022  "Expansion not of a proper type");
1023 
1024  int i, j;
1025  int cnt = 0;
1026  int nummodes0, nummodes1;
1027  int startvalue;
1028  if (outarray.size() != GetNcoeffs() - NumBndryCoeffs())
1029  {
1030  outarray = Array<OneD, unsigned int>(GetNcoeffs() - NumBndryCoeffs());
1031  }
1032 
1033  nummodes0 = m_base[0]->GetNumModes();
1034  nummodes1 = m_base[1]->GetNumModes();
1035 
1036  startvalue = 2 * nummodes1;
1037 
1038  for (i = 0; i < nummodes0 - 2; i++)
1039  {
1040  for (j = 0; j < nummodes1 - 3 - i; j++)
1041  {
1042  outarray[cnt++] = startvalue + j;
1043  }
1044  startvalue += nummodes1 - 2 - i;
1045  }
1046 }
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:130

References ASSERTL1, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), Nektar::StdRegions::StdExpansion::m_base, and Nektar::StdRegions::StdExpansion::NumBndryCoeffs().

Referenced by v_FwdTransBndConstrained().

◆ v_GetNtraces()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 747 of file StdTriExp.cpp.

748 {
749  return 3;
750 }

◆ v_GetNverts()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 742 of file StdTriExp.cpp.

743 {
744  return 3;
745 }

◆ v_GetSimplexEquiSpacedConnectivity()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1527 of file StdTriExp.cpp.

1529 {
1530  boost::ignore_unused(standard);
1531 
1532  int np1 = m_base[0]->GetNumPoints();
1533  int np2 = m_base[1]->GetNumPoints();
1534  int np = max(np1, np2);
1535 
1536  conn = Array<OneD, int>(3 * (np - 1) * (np - 1));
1537 
1538  int row = 0;
1539  int rowp1 = 0;
1540  int cnt = 0;
1541  for (int i = 0; i < np - 1; ++i)
1542  {
1543  rowp1 += np - i;
1544  for (int j = 0; j < np - i - 2; ++j)
1545  {
1546  conn[cnt++] = row + j;
1547  conn[cnt++] = row + j + 1;
1548  conn[cnt++] = rowp1 + j;
1549 
1550  conn[cnt++] = rowp1 + j + 1;
1551  conn[cnt++] = rowp1 + j;
1552  conn[cnt++] = row + j + 1;
1553  }
1554 
1555  conn[cnt++] = row + np - i - 2;
1556  conn[cnt++] = row + np - i - 1;
1557  conn[cnt++] = rowp1 + np - i - 2;
1558 
1559  row += np - i;
1560  }
1561 }

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

◆ v_GetTraceBasisKey()

const LibUtilities::BasisKey Nektar::StdRegions::StdTriExp::v_GetTraceBasisKey ( const int  i,
const int  j 
) const
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 857 of file StdTriExp.cpp.

859 {
860  boost::ignore_unused(j);
861  ASSERTL2(i >= 0 && i <= 2, "edge id is out of range");
862 
863  if (i == 0)
864  {
865  return GetBasis(0)->GetBasisKey();
866  }
867  else
868  {
869  switch (m_base[1]->GetBasisType())
870  {
872  {
873  switch (m_base[1]->GetPointsType())
874  {
875  case LibUtilities::eGaussRadauMAlpha1Beta0:
876  {
877  LibUtilities::PointsKey pkey(
878  m_base[1]
879  ->GetBasisKey()
880  .GetPointsKey()
881  .GetNumPoints() +
882  1,
884  return LibUtilities::BasisKey(LibUtilities::eModified_A,
885  m_base[1]->GetNumModes(),
886  pkey);
887  break;
888  }
889 
890  break;
891  // Currently this does not increase the points by
892  // 1 since when using this quadrature we are
893  // presuming it is already been increased by one
894  // when comopared to
895  // GaussRadauMAlpha1Beta0. Currently used in the
896  // GJP option
897  //
898  // Note have put down it back to numpoints +1 to
899  // test for use on tri faces and GJP.
900 
902  {
903  LibUtilities::PointsKey pkey(
904  m_base[1]
905  ->GetBasisKey()
906  .GetPointsKey()
907  .GetNumPoints() +
908  1,
910  return LibUtilities::BasisKey(LibUtilities::eModified_A,
911  m_base[1]->GetNumModes(),
912  pkey);
913  break;
914  }
915 
917  {
918  LibUtilities::PointsKey pkey(
919  m_base[1]
920  ->GetBasisKey()
921  .GetPointsKey()
922  .GetNumPoints() +
923  1,
925  return LibUtilities::BasisKey(LibUtilities::eModified_A,
926  m_base[1]->GetNumModes(),
927  pkey);
928  break;
929  }
930 
931  break;
932  default:
934  "Unexpected points distribution " +
936  [m_base[1]->GetPointsType()] +
937  " in StdTriExp::v_GetTraceBasisKey");
938  break;
939  }
940  break;
941  }
942  default:
944  "Information not available to set edge key");
945  break;
946  }
947  }
949 }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
const LibUtilities::BasisSharedPtr & GetBasis(int dir) const
This function gets the shared point to basis in the dir direction.
Definition: StdExpansion.h:118
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:211
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.
const std::string kPointsTypeStr[]
Definition: Foundations.hpp:54
@ eGaussRadauMLegendre
1D Gauss-Radau-Legendre quadrature points, pinned at x=-1
Definition: PointsType.h:49
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:53
@ ePolyEvenlySpaced
1D Evenly-spaced points using Lagrange polynomial
Definition: PointsType.h:75

References ASSERTL2, Nektar::ErrorUtil::efatal, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauMLegendre, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::LibUtilities::ePolyEvenlySpaced, Nektar::StdRegions::StdExpansion::GetBasis(), Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::GetNumPoints(), Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::LibUtilities::kPointsTypeStr, Nektar::StdRegions::StdExpansion::m_base, NEKERROR, and Nektar::LibUtilities::NullBasisKey().

◆ v_GetTraceCoeffMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion2D.

Definition at line 1080 of file StdTriExp.cpp.

1082 {
1083 
1086  "Mapping not defined for this type of basis");
1087 
1088  ASSERTL1(eid < 3, "eid must be between 0 and 2");
1089 
1090  int i;
1091  int order0 = m_base[0]->GetNumModes();
1092  int order1 = m_base[1]->GetNumModes();
1093  int numModes = (eid == 0) ? order0 : order1;
1094 
1095  if (maparray.size() != numModes)
1096  {
1097  maparray = Array<OneD, unsigned int>(numModes);
1098  }
1099 
1100  switch (eid)
1101  {
1102  case 0:
1103  {
1104  int cnt = 0;
1105  for (i = 0; i < numModes; cnt += order1 - i, ++i)
1106  {
1107  maparray[i] = cnt;
1108  }
1109  break;
1110  }
1111  case 1:
1112  {
1113  maparray[0] = order1;
1114  maparray[1] = 1;
1115  for (i = 2; i < numModes; i++)
1116  {
1117  maparray[i] = order1 - 1 + i;
1118  }
1119  break;
1120  }
1121  case 2:
1122  {
1123  for (i = 0; i < numModes; i++)
1124  {
1125  maparray[i] = i;
1126  }
1127  break;
1128  }
1129  default:
1130  ASSERTL0(false, "eid must be between 0 and 2");
1131  break;
1132  }
1133 }

References ASSERTL0, ASSERTL1, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::StdRegions::StdExpansion::GetBasisType(), and Nektar::StdRegions::StdExpansion::m_base.

◆ v_GetTraceInteriorToElementMap()

void Nektar::StdRegions::StdTriExp::v_GetTraceInteriorToElementMap ( const int  eid,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray,
const Orientation  edgeOrient = eForwards 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1135 of file StdTriExp.cpp.

1138 {
1141  "Mapping not defined for this type of basis");
1142  int i;
1143  const int nummodes1 = m_base[1]->GetNumModes();
1144  const int nEdgeIntCoeffs = GetTraceNcoeffs(eid) - 2;
1145 
1146  if (maparray.size() != nEdgeIntCoeffs)
1147  {
1148  maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
1149  }
1150 
1151  if (signarray.size() != nEdgeIntCoeffs)
1152  {
1153  signarray = Array<OneD, int>(nEdgeIntCoeffs, 1);
1154  }
1155  else
1156  {
1157  fill(signarray.get(), signarray.get() + nEdgeIntCoeffs, 1);
1158  }
1159 
1160  switch (eid)
1161  {
1162  case 0:
1163  {
1164  int cnt = 2 * nummodes1 - 1;
1165  for (i = 0; i < nEdgeIntCoeffs; cnt += nummodes1 - 2 - i, ++i)
1166  {
1167  maparray[i] = cnt;
1168  }
1169  break;
1170  }
1171  case 1:
1172  {
1173  for (i = 0; i < nEdgeIntCoeffs; i++)
1174  {
1175  maparray[i] = nummodes1 + 1 + i;
1176  }
1177  break;
1178  }
1179  case 2:
1180  {
1181  for (i = 0; i < nEdgeIntCoeffs; i++)
1182  {
1183  maparray[i] = 2 + i;
1184  }
1185  break;
1186  }
1187  default:
1188  {
1189  ASSERTL0(false, "eid must be between 0 and 2");
1190  break;
1191  }
1192  }
1193 
1194  if (edgeOrient == eBackwards)
1195  {
1196  for (i = 1; i < nEdgeIntCoeffs; i += 2)
1197  {
1198  signarray[i] = -1;
1199  }
1200  }
1201 }
int GetTraceNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th trace.
Definition: StdExpansion.h:267

References ASSERTL0, Nektar::StdRegions::eBackwards, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::GetTraceNcoeffs(), and Nektar::StdRegions::StdExpansion::m_base.

◆ v_GetTraceIntNcoeffs()

int Nektar::StdRegions::StdTriExp::v_GetTraceIntNcoeffs ( const int  i) const
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 791 of file StdTriExp.cpp.

792 {
793  ASSERTL2(i >= 0 && i <= 2, "edge id is out of range");
794 
795  if (i == 0)
796  {
797  return GetBasisNumModes(0) - 2;
798  }
799  else
800  {
801  return GetBasisNumModes(1) - 2;
802  }
803 }
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:175

References ASSERTL2, and Nektar::StdRegions::StdExpansion::GetBasisNumModes().

◆ v_GetTraceNcoeffs()

int Nektar::StdRegions::StdTriExp::v_GetTraceNcoeffs ( const int  i) const
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 777 of file StdTriExp.cpp.

778 {
779  ASSERTL2(i >= 0 && i <= 2, "edge id is out of range");
780 
781  if (i == 0)
782  {
783  return GetBasisNumModes(0);
784  }
785  else
786  {
787  return GetBasisNumModes(1);
788  }
789 }

References ASSERTL2, and Nektar::StdRegions::StdExpansion::GetBasisNumModes().

◆ v_GetTraceNumPoints()

int Nektar::StdRegions::StdTriExp::v_GetTraceNumPoints ( const int  i) const
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 805 of file StdTriExp.cpp.

806 {
807  ASSERTL2((i >= 0) && (i <= 2), "edge id is out of range");
808 
809  if (i == 0)
810  {
811  return GetNumPoints(0);
812  }
813  else
814  {
815  return GetNumPoints(1);
816  }
817 }

References ASSERTL2, and Nektar::StdRegions::StdExpansion::GetNumPoints().

◆ v_GetVertexMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp.

Definition at line 955 of file StdTriExp.cpp.

956 {
959  "Mapping not defined for this type of basis");
960 
961  int localDOF = 0;
962  if (useCoeffPacking == true)
963  {
964  switch (localVertexId)
965  {
966  case 0:
967  {
968  localDOF = 0;
969  break;
970  }
971  case 1:
972  {
973  localDOF = 1;
974  break;
975  }
976  case 2:
977  {
978  localDOF = m_base[1]->GetNumModes();
979  break;
980  }
981  default:
982  {
983  ASSERTL0(false, "eid must be between 0 and 2");
984  break;
985  }
986  }
987  }
988  else // follow book format for vertex indexing.
989  {
990  switch (localVertexId)
991  {
992  case 0:
993  {
994  localDOF = 0;
995  break;
996  }
997  case 1:
998  {
999  localDOF = m_base[1]->GetNumModes();
1000  break;
1001  }
1002  case 2:
1003  {
1004  localDOF = 1;
1005  break;
1006  }
1007  default:
1008  {
1009  ASSERTL0(false, "eid must be between 0 and 2");
1010  break;
1011  }
1012  }
1013  }
1014 
1015  return localDOF;
1016 }

References ASSERTL0, ASSERTL1, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::StdRegions::StdExpansion::GetBasisType(), and Nektar::StdRegions::StdExpansion::m_base.

◆ v_HelmholtzMatrixOp()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp, Nektar::LocalRegions::NodalTriExp, and Nektar::StdRegions::StdNodalTriExp.

Definition at line 1320 of file StdTriExp.cpp.

1323 {
1324  StdTriExp::v_HelmholtzMatrixOp_MatFree(inarray, outarray, mkey);
1325 }
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override

References Nektar::StdRegions::StdExpansion2D::v_HelmholtzMatrixOp_MatFree().

◆ v_Integral()

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

Integrates the specified function over the domain.

See also
StdRegions::StdExpansion::Integral.

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp, and Nektar::LocalRegions::NodalTriExp.

Definition at line 78 of file StdTriExp.cpp.

79 {
80  int i;
81  int nquad1 = m_base[1]->GetNumPoints();
82  Array<OneD, NekDouble> w1_tmp(nquad1);
83 
84  Array<OneD, const NekDouble> w0 = m_base[0]->GetW();
85  Array<OneD, const NekDouble> w1 = m_base[1]->GetW();
86  Array<OneD, const NekDouble> z1 = m_base[1]->GetZ();
87 
88  switch (m_base[1]->GetPointsType())
89  {
90  case LibUtilities::eGaussRadauMAlpha1Beta0: // (0,1) Jacobi Inner
91  // product
92  {
93  Vmath::Smul(nquad1, 0.5, w1, 1, w1_tmp, 1);
94  break;
95  }
96  default:
97  {
98  // include jacobian factor on whatever coordinates are defined.
99  for (i = 0; i < nquad1; ++i)
100  {
101  w1_tmp[i] = 0.5 * (1 - z1[i]) * w1[i];
102  }
103  break;
104  }
105  }
106 
107  return StdExpansion2D::Integral(inarray, w0, w1_tmp);
108 }
NekDouble Integral(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &w0, const Array< OneD, const NekDouble > &w1)

References Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::StdRegions::StdExpansion2D::Integral(), Nektar::StdRegions::StdExpansion::m_base, and Vmath::Smul().

◆ v_IProductWRTBase()

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

Calculate the inner product of inarray with respect to the basis B=base0[p]*base1[pq] and put into outarray.

\( \begin{array}{rcl} I_{pq} = (\phi^A_q \phi^B_{pq}, u) &=& \sum_{i=0}^{nq_0}\sum_{j=0}^{nq_1} \phi^A_p(\eta_{0,i})\phi^B_{pq}(\eta_{1,j}) w^0_i w^1_j u(\xi_{0,i} \xi_{1,j}) \\ & = & \sum_{i=0}^{nq_0} \phi^A_p(\eta_{0,i}) \sum_{j=0}^{nq_1} \phi^B_{pq}(\eta_{1,j}) \tilde{u}_{i,j} \end{array} \)

where

\( \tilde{u}_{i,j} = w^0_i w^1_j u(\xi_{0,i},\xi_{1,j}) \)

which can be implemented as

\( f_{pj} = \sum_{i=0}^{nq_0} \phi^A_p(\eta_{0,i}) \tilde{u}_{i,j} \rightarrow {\bf B_1 U} \) \( I_{pq} = \sum_{j=0}^{nq_1} \phi^B_{pq}(\eta_{1,j}) f_{pj} \rightarrow {\bf B_2[p*skip] f[skip]} \)

Recall: \( \eta_{1} = \frac{2(1+\xi_1)}{(1-\xi_2)}-1, \, \eta_2 = \xi_2\)

Note: For the orthgonality of this expansion to be realised the 'q' ordering must run fastest in contrast to the Quad and Hex ordering where 'p' index runs fastest to be consistent with the quadrature ordering.

In the triangular space the i (i.e. \(\eta_1\) direction) ordering still runs fastest by convention.

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp, and Nektar::LocalRegions::NodalTriExp.

Definition at line 444 of file StdTriExp.cpp.

446 {
447  StdTriExp::v_IProductWRTBase_SumFac(inarray, outarray);
448 }
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
Definition: StdTriExp.cpp:450

References v_IProductWRTBase_SumFac().

Referenced by v_FwdTrans(), and v_FwdTransBndConstrained().

◆ v_IProductWRTBase_SumFac()

void Nektar::StdRegions::StdTriExp::v_IProductWRTBase_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
bool  multiplybyweights = true 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp, Nektar::LocalRegions::TriExp, and Nektar::LocalRegions::NodalTriExp.

Definition at line 450 of file StdTriExp.cpp.

453 {
454  int nquad0 = m_base[0]->GetNumPoints();
455  int nquad1 = m_base[1]->GetNumPoints();
456  int order0 = m_base[0]->GetNumModes();
457 
458  if (multiplybyweights)
459  {
460  Array<OneD, NekDouble> tmp(nquad0 * nquad1 + nquad1 * order0);
461  Array<OneD, NekDouble> wsp(tmp + nquad0 * nquad1);
462 
463  // multiply by integration constants
464  MultiplyByQuadratureMetric(inarray, tmp);
465  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
466  m_base[1]->GetBdata(), tmp, outarray, wsp);
467  }
468  else
469  {
470  Array<OneD, NekDouble> wsp(nquad1 * order0);
471  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
472  m_base[1]->GetBdata(), inarray, outarray,
473  wsp);
474  }
475 }
void IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:729

References Nektar::StdRegions::StdExpansion2D::IProductWRTBase_SumFacKernel(), Nektar::StdRegions::StdExpansion::m_base, and Nektar::StdRegions::StdExpansion::MultiplyByQuadratureMetric().

Referenced by v_IProductWRTBase(), and Nektar::StdRegions::StdNodalTriExp::v_IProductWRTBase_SumFac().

◆ v_IProductWRTBase_SumFacKernel()

void Nektar::StdRegions::StdTriExp::v_IProductWRTBase_SumFacKernel ( const Array< OneD, const NekDouble > &  base0,
const Array< OneD, const NekDouble > &  base1,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
Array< OneD, NekDouble > &  wsp,
bool  doCheckCollDir0,
bool  doCheckCollDir1 
)
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion2D.

Definition at line 477 of file StdTriExp.cpp.

483 {
484  boost::ignore_unused(doCheckCollDir0, doCheckCollDir1);
485 
486  int i;
487  int mode;
488  int nquad0 = m_base[0]->GetNumPoints();
489  int nquad1 = m_base[1]->GetNumPoints();
490  int nmodes0 = m_base[0]->GetNumModes();
491  int nmodes1 = m_base[1]->GetNumModes();
492 
493  ASSERTL1(wsp.size() >= nquad1 * nmodes0,
494  "Workspace size is not sufficient");
495 
496  Blas::Dgemm('T', 'N', nquad1, nmodes0, nquad0, 1.0, inarray.get(), nquad0,
497  base0.get(), nquad0, 0.0, wsp.get(), nquad1);
498 
499  // Inner product with respect to 'b' direction
500  for (mode = i = 0; i < nmodes0; ++i)
501  {
502  Blas::Dgemv('T', nquad1, nmodes1 - i, 1.0, base1.get() + mode * nquad1,
503  nquad1, wsp.get() + i * nquad1, 1, 0.0,
504  outarray.get() + mode, 1);
505  mode += nmodes1 - i;
506  }
507 
508  // fix for modified basis by splitting top vertex mode
510  {
511  outarray[1] +=
512  Blas::Ddot(nquad1, base1.get() + nquad1, 1, wsp.get() + nquad1, 1);
513  }
514 }
static double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = .
Definition: Blas.hpp:182

References ASSERTL1, Blas::Ddot(), Blas::Dgemm(), Blas::Dgemv(), Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), and Nektar::StdRegions::StdExpansion::m_base.

◆ v_IProductWRTDerivBase()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp, Nektar::LocalRegions::TriExp, and Nektar::LocalRegions::NodalTriExp.

Definition at line 516 of file StdTriExp.cpp.

519 {
520  StdTriExp::v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
521 }
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdTriExp.cpp:523

References v_IProductWRTDerivBase_SumFac().

◆ v_IProductWRTDerivBase_SumFac()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp, Nektar::LocalRegions::TriExp, and Nektar::LocalRegions::NodalTriExp.

Definition at line 523 of file StdTriExp.cpp.

526 {
527  int i;
528  int nquad0 = m_base[0]->GetNumPoints();
529  int nquad1 = m_base[1]->GetNumPoints();
530  int nqtot = nquad0 * nquad1;
531  int nmodes0 = m_base[0]->GetNumModes();
532  int wspsize = max(max(nqtot, m_ncoeffs), nquad1 * nmodes0);
533 
534  Array<OneD, NekDouble> gfac0(2 * wspsize);
535  Array<OneD, NekDouble> tmp0(gfac0 + wspsize);
536 
537  const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
538 
539  // set up geometric factor: 2/(1-z1)
540  for (i = 0; i < nquad1; ++i)
541  {
542  gfac0[i] = 2.0 / (1 - z1[i]);
543  }
544 
545  for (i = 0; i < nquad1; ++i)
546  {
547  Vmath::Smul(nquad0, gfac0[i], &inarray[0] + i * nquad0, 1,
548  &tmp0[0] + i * nquad0, 1);
549  }
550 
551  MultiplyByQuadratureMetric(tmp0, tmp0);
552 
553  switch (dir)
554  {
555  case 0:
556  {
557  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
558  m_base[1]->GetBdata(), tmp0, outarray,
559  gfac0);
560  break;
561  }
562  case 1:
563  {
564  Array<OneD, NekDouble> tmp3(m_ncoeffs);
565  const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
566 
567  for (i = 0; i < nquad0; ++i)
568  {
569  gfac0[i] = 0.5 * (1 + z0[i]);
570  }
571 
572  for (i = 0; i < nquad1; ++i)
573  {
574  Vmath::Vmul(nquad0, &gfac0[0], 1, &tmp0[0] + i * nquad0, 1,
575  &tmp0[0] + i * nquad0, 1);
576  }
577 
578  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
579  m_base[1]->GetBdata(), tmp0, tmp3,
580  gfac0);
581 
582  MultiplyByQuadratureMetric(inarray, tmp0);
583  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
584  m_base[1]->GetDbdata(), tmp0, outarray,
585  gfac0);
586  Vmath::Vadd(m_ncoeffs, &tmp3[0], 1, &outarray[0], 1, &outarray[0],
587  1);
588  break;
589  }
590  default:
591  {
592  ASSERTL1(false, "input dir is out of range");
593  break;
594  }
595  }
596 }
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:359

References ASSERTL1, Nektar::StdRegions::StdExpansion2D::IProductWRTBase_SumFacKernel(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::StdRegions::StdExpansion::MultiplyByQuadratureMetric(), Vmath::Smul(), Vmath::Vadd(), and Vmath::Vmul().

Referenced by v_IProductWRTDerivBase(), and Nektar::StdRegions::StdNodalTriExp::v_IProductWRTDerivBase_SumFac().

◆ v_IsBoundaryInteriorExpansion()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 851 of file StdTriExp.cpp.

852 {
853  return m_base[0]->GetBasisType() == LibUtilities::eModified_A &&
854  m_base[1]->GetBasisType() == LibUtilities::eModified_B;
855 }

References Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, and Nektar::StdRegions::StdExpansion::m_base.

◆ v_LaplacianMatrixOp() [1/2]

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp, Nektar::LocalRegions::NodalTriExp, and Nektar::StdRegions::StdNodalTriExp.

Definition at line 1297 of file StdTriExp.cpp.

1300 {
1301  StdTriExp::v_LaplacianMatrixOp_MatFree(inarray, outarray, mkey);
1302 }
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override

References Nektar::StdRegions::StdExpansion2D::v_LaplacianMatrixOp_MatFree().

◆ v_LaplacianMatrixOp() [2/2]

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp, Nektar::LocalRegions::NodalTriExp, and Nektar::StdRegions::StdNodalTriExp.

Definition at line 1304 of file StdTriExp.cpp.

1308 {
1309  StdExpansion::LaplacianMatrixOp_MatFree(k1, k2, inarray, outarray, mkey);
1310 }
void LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)

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

◆ v_LocCollapsedToLocCoord()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 621 of file StdTriExp.cpp.

623 {
624  xi[0] = (1.0 + eta[0]) * (1.0 - eta[1]) * 0.5 - 1.0;
625  xi[1] = eta[1];
626 }

◆ v_LocCoordToLocCollapsed()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 602 of file StdTriExp.cpp.

604 {
605  NekDouble d1 = 1. - xi[1];
606  if (fabs(d1) < NekConstants::kNekZeroTol)
607  {
608  if (d1 >= 0.)
609  {
611  }
612  else
613  {
615  }
616  }
617  eta[0] = 2. * (1. + xi[0]) / d1 - 1.0;
618  eta[1] = xi[1];
619 }
static const NekDouble kNekZeroTol

References Nektar::NekConstants::kNekZeroTol.

◆ v_MassMatrixOp()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp, and Nektar::LocalRegions::NodalTriExp.

Definition at line 1290 of file StdTriExp.cpp.

1293 {
1294  StdExpansion::MassMatrixOp_MatFree(inarray, outarray, mkey);
1295 }
void MassMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)

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

◆ v_MultiplyByStdQuadratureMetric()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1487 of file StdTriExp.cpp.

1490 {
1491  int i;
1492  int nquad0 = m_base[0]->GetNumPoints();
1493  int nquad1 = m_base[1]->GetNumPoints();
1494 
1495  const Array<OneD, const NekDouble> &w0 = m_base[0]->GetW();
1496  const Array<OneD, const NekDouble> &w1 = m_base[1]->GetW();
1497  const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
1498 
1499  // multiply by integration constants
1500  for (i = 0; i < nquad1; ++i)
1501  {
1502  Vmath::Vmul(nquad0, inarray.get() + i * nquad0, 1, w0.get(), 1,
1503  outarray.get() + i * nquad0, 1);
1504  }
1505 
1506  switch (m_base[1]->GetPointsType())
1507  {
1508  // (1,0) Jacobi Inner product
1509  case LibUtilities::eGaussRadauMAlpha1Beta0:
1510  for (i = 0; i < nquad1; ++i)
1511  {
1512  Blas::Dscal(nquad0, 0.5 * w1[i], outarray.get() + i * nquad0,
1513  1);
1514  }
1515  break;
1516  // Legendre inner product
1517  default:
1518  for (i = 0; i < nquad1; ++i)
1519  {
1520  Blas::Dscal(nquad0, 0.5 * (1 - z1[i]) * w1[i],
1521  outarray.get() + i * nquad0, 1);
1522  }
1523  break;
1524  }
1525 }
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition: Blas.hpp:168

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

◆ v_NumBndryCoeffs()

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

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp.

Definition at line 757 of file StdTriExp.cpp.

758 {
760  "BasisType is not a boundary interior form");
762  "BasisType is not a boundary interior form");
763 
764  return 3 + (GetBasisNumModes(0) - 2) + 2 * (GetBasisNumModes(1) - 2);
765 }

References ASSERTL1, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::StdRegions::StdExpansion::GetBasisNumModes(), and Nektar::StdRegions::StdExpansion::GetBasisType().

Referenced by v_FwdTransBndConstrained().

◆ v_NumDGBndryCoeffs()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 767 of file StdTriExp.cpp.

768 {
770  "BasisType is not a boundary interior form");
772  "BasisType is not a boundary interior form");
773 
774  return GetBasisNumModes(0) + 2 * GetBasisNumModes(1);
775 }

References ASSERTL1, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::StdRegions::StdExpansion::GetBasisNumModes(), and Nektar::StdRegions::StdExpansion::GetBasisType().

◆ v_PhysDeriv() [1/2]

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

Calculate the derivative of the physical points.

\( \frac{\partial u}{\partial x_1} = \left . \frac{2.0}{1-\eta_2} \frac{\partial u}{\partial d\eta_1} \right |_{\eta_2}\)

\( \frac{\partial u}{\partial x_2} = \left . \frac{1+\eta_1}{1-\eta_2} \frac{\partial u}{\partial d\eta_1} \right |_{\eta_2} + \left . \frac{\partial u}{\partial d\eta_2} \right |_{\eta_1} \)

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp, and Nektar::LocalRegions::NodalTriExp.

Definition at line 126 of file StdTriExp.cpp.

130 {
131  boost::ignore_unused(out_d2);
132 
133  int i;
134  int nquad0 = m_base[0]->GetNumPoints();
135  int nquad1 = m_base[1]->GetNumPoints();
136  Array<OneD, NekDouble> wsp(std::max(nquad0, nquad1));
137 
138  const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
139  const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
140 
141  // set up geometric factor: 2/(1-z1)
142  Vmath::Sadd(nquad1, -1.0, z1, 1, wsp, 1);
143  Vmath::Sdiv(nquad1, -2.0, wsp, 1, wsp, 1);
144 
145  if (out_d0.size() > 0)
146  {
147  PhysTensorDeriv(inarray, out_d0, out_d1);
148 
149  for (i = 0; i < nquad1; ++i)
150  {
151  Blas::Dscal(nquad0, wsp[i], &out_d0[0] + i * nquad0, 1);
152  }
153 
154  // if no d1 required do not need to calculate both deriv
155  if (out_d1.size() > 0)
156  {
157  // set up geometric factor: (1+z0)/2
158  Vmath::Sadd(nquad0, 1.0, z0, 1, wsp, 1);
159  Vmath::Smul(nquad0, 0.5, wsp, 1, wsp, 1);
160 
161  for (i = 0; i < nquad1; ++i)
162  {
163  Vmath::Vvtvp(nquad0, &wsp[0], 1, &out_d0[0] + i * nquad0, 1,
164  &out_d1[0] + i * nquad0, 1,
165  &out_d1[0] + i * nquad0, 1);
166  }
167  }
168  }
169  else if (out_d1.size() > 0)
170  {
171  Array<OneD, NekDouble> diff0(nquad0 * nquad1);
172  PhysTensorDeriv(inarray, diff0, out_d1);
173 
174  for (i = 0; i < nquad1; ++i)
175  {
176  Blas::Dscal(nquad0, wsp[i], &diff0[0] + i * nquad0, 1);
177  }
178 
179  Vmath::Sadd(nquad0, 1.0, z0, 1, wsp, 1);
180  Vmath::Smul(nquad0, 0.5, wsp, 1, wsp, 1);
181 
182  for (i = 0; i < nquad1; ++i)
183  {
184  Vmath::Vvtvp(nquad0, &wsp[0], 1, &diff0[0] + i * nquad0, 1,
185  &out_d1[0] + i * nquad0, 1, &out_d1[0] + i * nquad0,
186  1);
187  }
188  }
189 }
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d0, Array< OneD, NekDouble > &outarray_d1)
Calculate the 2D derivative in the local tensor/collapsed coordinate at the physical points.
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:574
void Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/y.
Definition: Vmath.cpp:324
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add scalar y = alpha + x.
Definition: Vmath.cpp:384

References Blas::Dscal(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion2D::PhysTensorDeriv(), Vmath::Sadd(), Vmath::Sdiv(), Vmath::Smul(), and Vmath::Vvtvp().

Referenced by v_PhysDeriv(), and v_StdPhysDeriv().

◆ v_PhysDeriv() [2/2]

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

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

See also
StdRegions::StdExpansion::PhysDeriv

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp, and Nektar::LocalRegions::NodalTriExp.

Definition at line 191 of file StdTriExp.cpp.

194 {
195  switch (dir)
196  {
197  case 0:
198  {
199  v_PhysDeriv(inarray, outarray, NullNekDouble1DArray);
200  break;
201  }
202  case 1:
203  {
204  v_PhysDeriv(inarray, NullNekDouble1DArray, outarray);
205  break;
206  }
207  default:
208  {
209  ASSERTL1(false, "input dir is out of range");
210  break;
211  }
212  }
213 }
void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray) override
Calculate the derivative of the physical points.
Definition: StdTriExp.cpp:126
static Array< OneD, NekDouble > NullNekDouble1DArray

References ASSERTL1, Nektar::NullNekDouble1DArray, and v_PhysDeriv().

◆ v_PhysEvaluate()

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

Reimplemented from Nektar::StdRegions::StdExpansion2D.

Reimplemented in Nektar::LocalRegions::TriExp, and Nektar::LocalRegions::NodalTriExp.

Definition at line 697 of file StdTriExp.cpp.

700 {
701  // Collapse coordinates
702  Array<OneD, NekDouble> coll(2, 0.0);
703  LocCoordToLocCollapsed(coord, coll);
704 
705  // If near singularity do the old interpolation matrix method
706  if ((1 - coll[1]) < 1e-5)
707  {
708  int totPoints = GetTotPoints();
709  Array<OneD, NekDouble> EphysDeriv0(totPoints), EphysDeriv1(totPoints);
710  PhysDeriv(inarray, EphysDeriv0, EphysDeriv1);
711 
712  Array<OneD, DNekMatSharedPtr> I(2);
713  I[0] = GetBase()[0]->GetI(coll);
714  I[1] = GetBase()[1]->GetI(coll + 1);
715 
716  firstOrderDerivs[0] = PhysEvaluate(I, EphysDeriv0);
717  firstOrderDerivs[1] = PhysEvaluate(I, EphysDeriv1);
718  return PhysEvaluate(I, inarray);
719  }
720 
721  // set up geometric factor: 2.0/(1.0-z1)
722  NekDouble fac0 = 2 / (1 - coll[1]);
723 
724  NekDouble val = BaryTensorDeriv(coll, inarray, firstOrderDerivs);
725 
726  // Copy d0 into temp for d1
727  NekDouble temp;
728  temp = firstOrderDerivs[0];
729 
730  // Multiply by geometric factor
731  firstOrderDerivs[0] = firstOrderDerivs[0] * fac0;
732 
733  // set up geometric factor: (1+z0)/(1-z1)
734  NekDouble fac1 = fac0 * (coll[0] + 1) / 2;
735 
736  // Multiply out_d0 by geometric factor and add to out_d1
737  firstOrderDerivs[1] += fac1 * temp;
738 
739  return val;
740 }
NekDouble BaryTensorDeriv(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140
const Array< OneD, const LibUtilities::BasisSharedPtr > & GetBase() const
This function gets the shared point to basis.
Definition: StdExpansion.h:106
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.
Definition: StdExpansion.h:918
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)
Definition: StdExpansion.h:848

References Nektar::StdRegions::StdExpansion2D::BaryTensorDeriv(), Nektar::StdRegions::StdExpansion::GetBase(), Nektar::StdRegions::StdExpansion::GetTotPoints(), Nektar::StdRegions::StdExpansion::LocCoordToLocCollapsed(), Nektar::StdRegions::StdExpansion::PhysDeriv(), and Nektar::StdRegions::StdExpansion::PhysEvaluate().

◆ v_PhysEvaluateBasis()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 673 of file StdTriExp.cpp.

675 {
676  Array<OneD, NekDouble> coll(2);
677  LocCoordToLocCollapsed(coords, coll);
678 
679  // From mode we need to determine mode0 and mode1 in the (p,q)
680  // direction. mode1 can be directly inferred from mode.
681  const int nm1 = m_base[1]->GetNumModes();
682  const double c = 1 + 2 * nm1;
683  const int mode0 = floor(0.5 * (c - sqrt(c * c - 8 * mode)));
684 
685  if (mode == 1 && m_base[0]->GetBasisType() == LibUtilities::eModified_A)
686  {
687  // Account for collapsed vertex.
688  return StdExpansion::BaryEvaluateBasis<1>(coll[1], 1);
689  }
690  else
691  {
692  return StdExpansion::BaryEvaluateBasis<0>(coll[0], mode0) *
693  StdExpansion::BaryEvaluateBasis<1>(coll[1], mode);
694  }
695 }
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::LocCoordToLocCollapsed(), Nektar::StdRegions::StdExpansion::m_base, and tinysimd::sqrt().

◆ v_ReduceOrderCoeffs()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp.

Definition at line 1429 of file StdTriExp.cpp.

1432 {
1433  int n_coeffs = inarray.size();
1434  int nquad0 = m_base[0]->GetNumPoints();
1435  int nquad1 = m_base[1]->GetNumPoints();
1436  Array<OneD, NekDouble> coeff(n_coeffs);
1437  Array<OneD, NekDouble> coeff_tmp(n_coeffs, 0.0);
1438  Array<OneD, NekDouble> tmp;
1439  Array<OneD, NekDouble> tmp2;
1440  int nqtot = nquad0 * nquad1;
1441  Array<OneD, NekDouble> phys_tmp(nqtot, 0.0);
1442 
1443  int nmodes0 = m_base[0]->GetNumModes();
1444  int nmodes1 = m_base[1]->GetNumModes();
1445  int numMin2 = nmodes0;
1446  int i;
1447 
1448  const LibUtilities::PointsKey Pkey0(nmodes0,
1450  const LibUtilities::PointsKey Pkey1(nmodes1,
1452 
1453  LibUtilities::BasisKey b0(m_base[0]->GetBasisType(), nmodes0, Pkey0);
1454  LibUtilities::BasisKey b1(m_base[1]->GetBasisType(), nmodes1, Pkey1);
1455 
1456  LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A, nmodes0, Pkey0);
1457  LibUtilities::BasisKey bortho1(LibUtilities::eOrtho_B, nmodes1, Pkey1);
1458 
1459  StdRegions::StdTriExpSharedPtr m_OrthoTriExp;
1461 
1464  bortho0, bortho1);
1465 
1466  m_TriExp->BwdTrans(inarray, phys_tmp);
1467  m_OrthoTriExp->FwdTrans(phys_tmp, coeff);
1468 
1469  for (i = 0; i < n_coeffs; i++)
1470  {
1471  if (i == numMin)
1472  {
1473  coeff[i] = 0.0;
1474  numMin += numMin2 - 1;
1475  numMin2 -= 1.0;
1476  }
1477  }
1478 
1479  m_OrthoTriExp->BwdTrans(coeff, phys_tmp);
1480  m_TriExp->FwdTrans(phys_tmp, outarray);
1481 }
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:44
std::shared_ptr< StdTriExp > StdTriExpSharedPtr
Definition: StdTriExp.h:232

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::eOrtho_B, Nektar::StdRegions::StdExpansion::GetBasisType(), and Nektar::StdRegions::StdExpansion::m_base.

◆ v_StdPhysDeriv() [1/2]

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::NodalTriExp.

Definition at line 215 of file StdTriExp.cpp.

219 {
220  boost::ignore_unused(out_d2);
221  StdTriExp::v_PhysDeriv(inarray, out_d0, out_d1);
222 }

References v_PhysDeriv().

◆ v_StdPhysDeriv() [2/2]

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 224 of file StdTriExp.cpp.

227 {
228  StdTriExp::v_PhysDeriv(dir, inarray, outarray);
229 }

References v_PhysDeriv().

◆ v_SVVLaplacianFilter()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp.

Definition at line 1327 of file StdTriExp.cpp.

1329 {
1330  int qa = m_base[0]->GetNumPoints();
1331  int qb = m_base[1]->GetNumPoints();
1332  int nmodes_a = m_base[0]->GetNumModes();
1333  int nmodes_b = m_base[1]->GetNumModes();
1334 
1335  // Declare orthogonal basis.
1336  LibUtilities::PointsKey pa(qa, m_base[0]->GetPointsType());
1337  LibUtilities::PointsKey pb(qb, m_base[1]->GetPointsType());
1338 
1339  LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A, nmodes_a, pa);
1340  LibUtilities::BasisKey Bb(LibUtilities::eOrtho_B, nmodes_b, pb);
1341  StdTriExp OrthoExp(Ba, Bb);
1342 
1343  Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
1344 
1345  // project onto physical space.
1346  OrthoExp.FwdTrans(array, orthocoeffs);
1347 
1348  if (mkey.ConstFactorExists(
1349  eFactorSVVPowerKerDiffCoeff)) // Rodrigo's power kern
1350  {
1351  NekDouble cutoff = mkey.GetConstFactor(eFactorSVVCutoffRatio);
1352  NekDouble SvvDiffCoeff =
1353  mkey.GetConstFactor(eFactorSVVPowerKerDiffCoeff) *
1354  mkey.GetConstFactor(eFactorSVVDiffCoeff);
1355 
1356  int cnt = 0;
1357  for (int j = 0; j < nmodes_a; ++j)
1358  {
1359  for (int k = 0; k < nmodes_b - j; ++k, ++cnt)
1360  {
1361  NekDouble fac = std::max(
1362  pow((1.0 * j) / (nmodes_a - 1), cutoff * nmodes_a),
1363  pow((1.0 * k) / (nmodes_b - 1), cutoff * nmodes_b));
1364 
1365  orthocoeffs[cnt] *= (SvvDiffCoeff * fac);
1366  }
1367  }
1368  }
1369  else if (mkey.ConstFactorExists(
1370  eFactorSVVDGKerDiffCoeff)) // Rodrigo/mansoor's DG kernel
1371  {
1372  NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDGKerDiffCoeff) *
1373  mkey.GetConstFactor(eFactorSVVDiffCoeff);
1374  int max_ab = max(nmodes_a - kSVVDGFiltermodesmin,
1375  nmodes_b - kSVVDGFiltermodesmin);
1376  max_ab = max(max_ab, 0);
1377  max_ab = min(max_ab, kSVVDGFiltermodesmax - kSVVDGFiltermodesmin);
1378 
1379  int cnt = 0;
1380  for (int j = 0; j < nmodes_a; ++j)
1381  {
1382  for (int k = 0; k < nmodes_b - j; ++k, ++cnt)
1383  {
1384  int maxjk = max(j, k);
1385  maxjk = min(maxjk, kSVVDGFiltermodesmax - 1);
1386 
1387  orthocoeffs[cnt] *= SvvDiffCoeff * kSVVDGFilter[max_ab][maxjk];
1388  }
1389  }
1390  }
1391  else
1392  {
1393  NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDiffCoeff);
1394 
1395  int cutoff = (int)(mkey.GetConstFactor(eFactorSVVCutoffRatio) *
1396  min(nmodes_a, nmodes_b));
1397 
1398  NekDouble epsilon = 1.0;
1399  int nmodes = min(nmodes_a, nmodes_b);
1400 
1401  int cnt = 0;
1402 
1403  // apply SVV filter (JEL)
1404  for (int j = 0; j < nmodes_a; ++j)
1405  {
1406  for (int k = 0; k < nmodes_b - j; ++k)
1407  {
1408  if (j + k >= cutoff)
1409  {
1410  orthocoeffs[cnt] *=
1411  (SvvDiffCoeff *
1412  exp(-(j + k - nmodes) * (j + k - nmodes) /
1413  ((NekDouble)((j + k - cutoff + epsilon) *
1414  (j + k - cutoff + epsilon)))));
1415  }
1416  else
1417  {
1418  orthocoeffs[cnt] *= 0.0;
1419  }
1420  cnt++;
1421  }
1422  }
1423  }
1424 
1425  // backward transform to physical space
1426  OrthoExp.BwdTrans(orthocoeffs, array);
1427 }
const int kSVVDGFiltermodesmin
Definition: StdRegions.hpp:469
const int kSVVDGFiltermodesmax
Definition: StdRegions.hpp:470
const NekDouble kSVVDGFilter[9][11]
Definition: StdRegions.hpp:472

References Nektar::StdRegions::StdExpansion::BwdTrans(), Nektar::StdRegions::StdMatrixKey::ConstFactorExists(), Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::eFactorSVVDGKerDiffCoeff, Nektar::StdRegions::eFactorSVVDiffCoeff, Nektar::StdRegions::eFactorSVVPowerKerDiffCoeff, Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::eOrtho_B, Nektar::StdRegions::StdExpansion::FwdTrans(), Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::StdRegions::kSVVDGFilter, Nektar::StdRegions::kSVVDGFiltermodesmax, Nektar::StdRegions::kSVVDGFiltermodesmin, and Nektar::StdRegions::StdExpansion::m_base.

◆ v_WeakDerivMatrixOp()

void Nektar::StdRegions::StdTriExp::v_WeakDerivMatrixOp ( const int  i,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdMatrixKey mkey 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp, Nektar::LocalRegions::NodalTriExp, and Nektar::StdRegions::StdNodalTriExp.

Definition at line 1312 of file StdTriExp.cpp.

1316 {
1317  StdExpansion::WeakDerivMatrixOp_MatFree(i, inarray, outarray, mkey);
1318 }
void WeakDerivMatrixOp_MatFree(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)

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