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)
 
 ~StdTriExp ()
 
- 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 ()
 
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)
 
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 (void) const
 
std::shared_ptr< StdExpansionGetLinStdExp (void) const
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion ()
 
bool IsNodalNonTensorialExp ()
 
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space. More...
 
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space. More...
 
void FwdTrans_BndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble Integral (const Array< OneD, const NekDouble > &inarray)
 This function integrates the specified function over the domain. More...
 
void FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 This function fills the array outarray with the mode-th mode of the expansion. More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 this function calculates the inner product of a given function f with the different modes of the expansion More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check)
 
void IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
int GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void SetElmtId (const int id)
 Set the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2=NullNekDouble1DArray, Array< OneD, NekDouble > &coords_3=NullNekDouble1DArray)
 this function returns the physical coordinates of the quadrature points of the expansion More...
 
void GetCoord (const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
 given the coordinates of a point of the element in the local collapsed coordinate system, this function calculates the physical coordinates of the point More...
 
DNekMatSharedPtr GetStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr GetStdStaticCondMatrix (const StdMatrixKey &mkey)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
void DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
int CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
NekDouble StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
int GetCoordim ()
 
void GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
void GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
int GetVertexMap (const int localVertexId, bool useCoeffPacking=false)
 
void GetTraceToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
 
void GetTraceInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eForwards)
 
void GetTraceNumModes (const int tid, int &numModes0, int &numModes1, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2)
 
void MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
DNekMatSharedPtr CreateGeneralMatrix (const StdMatrixKey &mkey)
 this function generates the mass matrix \(\mathbf{M}[i][j] = \int \phi_i(\mathbf{x}) \phi_j(\mathbf{x}) d\mathbf{x}\) More...
 
void GeneralMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
void ExponentialFilter (Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
 
void LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
DNekMatSharedPtr GenMatrix (const StdMatrixKey &mkey)
 
void PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void PhysDeriv_s (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds)
 
void PhysDeriv_n (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn)
 
void PhysDirectionalDeriv (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &outarray)
 
void StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluateBasis (const Array< OneD, const NekDouble > &coords, int mode)
 This function evaluates the basis function mode mode at a point coords of the domain. More...
 
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta. More...
 
void LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi)
 Convert local collapsed coordinates eta into local cartesian coordinate xi. More...
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual void v_DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
NekDouble Linf (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( L_\infty\) error \( |\epsilon|_\infty = \max |u - u_{exact}|\) where \( u_{exact}\) is given by the array sol. More...
 
NekDouble L2 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( L_2\) error, \( | \epsilon |_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 dx \right]^{1/2} d\xi_1 \) where \( u_{exact}\) is given by the array sol. More...
 
NekDouble H1 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( H^1\) error, \( | \epsilon |^1_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 + \nabla(u - u_{exact})\cdot\nabla(u - u_{exact})\cdot dx \right]^{1/2} d\xi_1 \) where \( u_{exact}\) is given by the array sol. More...
 
const LibUtilities::PointsKeyVector GetPointsKeys () const
 
DNekMatSharedPtr BuildInverseTransformationMatrix (const DNekScalMatSharedPtr &m_transformationmatrix)
 
void PhysInterpToSimplexEquiSpaced (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int npset=-1)
 This function performs an interpolation from the physical space points provided at input into an array of equispaced points which are not the collapsed coordinate. So for a tetrahedron you will only get a tetrahedral number of values. More...
 
void GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true)
 This function provides the connectivity of local simplices (triangles or tets) to connect the equispaced data points provided by PhysInterpToSimplexEquiSpaced. More...
 
void EquiSpacedToCoeffs (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs a projection/interpolation from the equispaced points sometimes used in post-processing onto the coefficient space. More...
 
template<class T >
std::shared_ptr< T > as ()
 
void IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
void GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 

Protected Member Functions

virtual NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray)
 Integrates the specified function over the domain. More...
 
virtual 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)
 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)
 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)
 
virtual void v_StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Backward tranform for triangular elements. More...
 
virtual void v_BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
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)
 
virtual void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Transform a given function from physical quadrature space to coefficient space. More...
 
virtual void v_FwdTrans_BndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 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_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
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)
 
virtual void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTDerivBase_MatOp (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 
virtual void v_LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi)
 
virtual void v_FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 
NekDouble v_PhysEvaluateBasis (const Array< OneD, const NekDouble > &coords, int mode) final
 
virtual int v_GetNverts () const
 
virtual int v_GetNtraces () const
 
virtual LibUtilities::ShapeType v_DetShapeType () const
 
virtual int v_NumBndryCoeffs () const
 
virtual int v_NumDGBndryCoeffs () const
 
virtual int v_GetTraceNcoeffs (const int i) const
 
virtual int v_GetTraceNumPoints (const int i) const
 
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
virtual void v_GetCoords (Array< OneD, NekDouble > &coords_x, Array< OneD, NekDouble > &coords_y, Array< OneD, NekDouble > &coords_z)
 
virtual bool v_IsBoundaryInteriorExpansion ()
 
virtual const LibUtilities::BasisKey v_GetTraceBasisKey (const int i, const int j) const
 
virtual int v_GetVertexMap (int localVertexId, bool useCoeffPacking=false)
 
virtual void v_GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
virtual void v_GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
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)
 
virtual void v_GetTraceInteriorToElementMap (const int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation edgeOrient=eForwards)
 
virtual DNekMatSharedPtr v_GenMatrix (const StdMatrixKey &mkey)
 
virtual DNekMatSharedPtr v_CreateStdMatrix (const StdMatrixKey &mkey)
 
virtual void v_MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
virtual void v_ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_GeneralMatrixOp_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true)
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion2D
virtual NekDouble v_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...
 
virtual NekDouble v_PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 
virtual void v_LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
 
virtual void v_HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
 
virtual void v_GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion
DNekMatSharedPtr CreateStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr CreateStdStaticCondMatrix (const StdMatrixKey &mkey)
 Create the static condensation of a matrix when using a boundary interior decomposition. More...
 
void BwdTrans_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase_SumFac (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GeneralMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
 
void LaplacianMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp_MatFree (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void HelmholtzMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
 
virtual NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
template<int DIR>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals)
 This function performs the barycentric interpolation of the polynomial stored in coord at a point physvals using barycentric interpolation weights in direction. More...
 
template<int DIR>
NekDouble BaryEvaluateBasis (const NekDouble &coord, const int &mode)
 

Additional Inherited Members

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

Detailed Description

Definition at line 49 of file StdTriExp.h.

Constructor & Destructor Documentation

◆ StdTriExp() [1/3]

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

Definition at line 49 of file StdTriExp.cpp.

50  {
51  }

◆ StdTriExp() [2/3]

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

Definition at line 54 of file StdTriExp.cpp.

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

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

◆ StdTriExp() [3/3]

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

Definition at line 71 of file StdTriExp.cpp.

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

◆ ~StdTriExp()

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

Definition at line 77 of file StdTriExp.cpp.

78  {
79  }

Member Function Documentation

◆ v_BwdTrans()

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

Backward tranform for triangular elements.

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

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp.

Definition at line 251 of file StdTriExp.cpp.

254  {
255  v_BwdTrans_SumFac(inarray,outarray);
256  }
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdTriExp.cpp:259

References v_BwdTrans_SumFac().

◆ v_BwdTrans_SumFac()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 259 of file StdTriExp.cpp.

262  {
263  Array<OneD, NekDouble> wsp(m_base[0]->GetNumPoints()*
264  m_base[1]->GetNumModes());
265 
266  BwdTrans_SumFacKernel(m_base[0]->GetBdata(),
267  m_base[1]->GetBdata(),
268  inarray,outarray,wsp);
269  }
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:221
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 
)
protectedvirtual

Implements Nektar::StdRegions::StdExpansion2D.

Definition at line 271 of file StdTriExp.cpp.

279  {
280  boost::ignore_unused(doCheckCollDir0, doCheckCollDir1);
281 
282  int i;
283  int mode;
284  int nquad0 = m_base[0]->GetNumPoints();
285  int nquad1 = m_base[1]->GetNumPoints();
286  int nmodes0 = m_base[0]->GetNumModes();
287  int nmodes1 = m_base[1]->GetNumModes();
288 
289  ASSERTL1(wsp.size() >= nquad0*nmodes1,
290  "Workspace size is not sufficient");
293  "Basis[1] is not of general tensor type");
294 
295  for (i = mode = 0; i < nmodes0; ++i)
296  {
297  Blas::Dgemv('N', nquad1,nmodes1-i,1.0,base1.get()+mode*nquad1,
298  nquad1,&inarray[0]+mode,1,0.0,&wsp[0]+i*nquad1,1);
299  mode += nmodes1-i;
300  }
301 
302  // fix for modified basis by splitting top vertex mode
304  {
305  Blas::Daxpy(nquad1,inarray[1],base1.get()+nquad1,1,
306  &wsp[0]+nquad1,1);
307  }
308 
309  Blas::Dgemm('N','T', nquad0,nquad1,nmodes0,1.0, base0.get(),nquad0,
310  &wsp[0], nquad1,0.0, &outarray[0], nquad0);
311  }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:274
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:158
static void Dgemv(const char &trans, const int &m, const int &n, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = A x where A[m x n].
Definition: Blas.hpp:265
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:394
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
Definition: Blas.hpp:167
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:49
@ eOrtho_B
Principle Orthogonal Functions .
Definition: BasisType.h:46
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 841 of file StdTriExp.cpp.

844  {
846  nummodes[modes_offset],
847  nummodes[modes_offset+1]);
848  modes_offset += 2;
849 
850  return nmodes;
851  }

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

◆ v_CreateStdMatrix()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 1323 of file StdTriExp.cpp.

1324  {
1325  return v_GenMatrix(mkey);
1326  }
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey)
Definition: StdTriExp.cpp:1246

References v_GenMatrix().

◆ v_DetShapeType()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 788 of file StdTriExp.cpp.

789  {
791  }

References Nektar::LibUtilities::eTriangle.

◆ v_FillMode()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp.

Definition at line 703 of file StdTriExp.cpp.

705  {
706  int i,m;
707  int nquad0 = m_base[0]->GetNumPoints();
708  int nquad1 = m_base[1]->GetNumPoints();
709  int order0 = m_base[0]->GetNumModes();
710  int order1 = m_base[1]->GetNumModes();
711  int mode0 = 0;
712  Array<OneD, const NekDouble> base0 = m_base[0]->GetBdata();
713  Array<OneD, const NekDouble> base1 = m_base[1]->GetBdata();
714 
715  ASSERTL2(mode <= m_ncoeffs,
716  "calling argument mode is larger than "
717  "total expansion order");
718 
719  m = order1;
720  for (i = 0; i < order0; ++i, m+=order1-i)
721  {
722  if (m > mode)
723  {
724  mode0 = i;
725  break;
726  }
727  }
728 
729  // deal with top vertex mode in modified basis
730  if (mode == 1 &&
732  {
733  Vmath::Fill(nquad0*nquad1 , 1.0, outarray, 1);
734  }
735  else
736  {
737  for (i = 0; i < nquad1; ++i)
738  {
739  Vmath::Vcopy(nquad0,(NekDouble *)(base0.get()+mode0*nquad0),
740  1,&outarray[0]+i*nquad0,1);
741  }
742  }
743 
744  for(i = 0; i < nquad0; ++i)
745  {
746  Vmath::Vmul(nquad1,(NekDouble *)(base1.get() + mode*nquad1),
747  1,&outarray[0]+i,nquad0,&outarray[0]+i,nquad0);
748  }
749  }
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:192
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:1199

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

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 313 of file StdTriExp.cpp.

316  {
317  v_IProductWRTBase(inarray,outarray);
318 
319  // get Mass matrix inverse
320  StdMatrixKey masskey(eInvMass,DetShapeType(),*this);
321  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
322 
323  // copy inarray in case inarray == outarray
324  NekVector<NekDouble> in(m_ncoeffs,outarray,eCopy);
325  NekVector<NekDouble> out(m_ncoeffs,outarray,eWrapper);
326 
327  out = (*matsys)*in;
328  }
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:617
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:376
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Calculate the inner product of inarray with respect to the basis B=base0[p]*base1[pq] and put into ou...
Definition: StdTriExp.cpp:466
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp.

Definition at line 331 of file StdTriExp.cpp.

334  {
335  int i,j;
336  int npoints[2] = {m_base[0]->GetNumPoints(),
337  m_base[1]->GetNumPoints()};
338  int nmodes[2] = {m_base[0]->GetNumModes(),
339  m_base[1]->GetNumModes()};
340 
341  fill(outarray.get(), outarray.get()+m_ncoeffs, 0.0 );
342 
343  Array<OneD, NekDouble> physEdge[3];
344  Array<OneD, NekDouble> coeffEdge[3];
345  for(i = 0; i < 3; i++)
346  {
347  physEdge[i] = Array<OneD, NekDouble>(npoints[i!=0]);
348  coeffEdge[i] = Array<OneD, NekDouble>(nmodes[i!=0]);
349  }
350 
351  for(i = 0; i < npoints[0]; i++)
352  {
353  physEdge[0][i] = inarray[i];
354  }
355 
356  for(i = 0; i < npoints[1]; i++)
357  {
358  physEdge[1][i] = inarray[npoints[0]-1+i*npoints[0]];
359  physEdge[2][i] = inarray[(npoints[1]-1)*npoints[0]-i*npoints[0]];
360  }
361 
362  StdSegExpSharedPtr segexp[2] = {
364  m_base[0]->GetBasisKey()),
366  m_base[1]->GetBasisKey())
367  };
368 
369  Array<OneD, unsigned int> mapArray;
370  Array<OneD, int> signArray;
371  NekDouble sign;
372 
373  for (i = 0; i < 3; i++)
374  {
375  segexp[i!=0]->FwdTrans_BndConstrained(physEdge[i],coeffEdge[i]);
376 
377  GetTraceToElementMap(i,mapArray,signArray);
378  for (j = 0; j < nmodes[i != 0]; j++)
379  {
380  sign = (NekDouble) signArray[j];
381  outarray[ mapArray[j] ] = sign * coeffEdge[i][j];
382  }
383  }
384 
385  Array<OneD, NekDouble> tmp0(m_ncoeffs);
386  Array<OneD, NekDouble> tmp1(m_ncoeffs);
387 
388  StdMatrixKey masskey(eMass,DetShapeType(),*this);
389  MassMatrixOp(outarray,tmp0,masskey);
390  v_IProductWRTBase(inarray,tmp1);
391 
392  Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
393 
394  // get Mass matrix inverse (only of interior DOF)
395  // use block (1,1) of the static condensed system
396  // note: this block alreay contains the inverse matrix
397  DNekMatSharedPtr matsys =
398  (m_stdStaticCondMatrixManager[masskey])->GetBlock(1,1);
399 
400  int nBoundaryDofs = v_NumBndryCoeffs();
401  int nInteriorDofs = m_ncoeffs - nBoundaryDofs;
402 
403  Array<OneD, NekDouble> rhs (nInteriorDofs);
404  Array<OneD, NekDouble> result(nInteriorDofs);
405 
406  v_GetInteriorMap(mapArray);
407 
408  for (i = 0; i < nInteriorDofs; i++)
409  {
410  rhs[i] = tmp1[ mapArray[i] ];
411  }
412 
413  Blas::Dgemv('N',nInteriorDofs,nInteriorDofs,
414  1.0,&(matsys->GetPtr())[0],nInteriorDofs,
415  rhs.get(),1,
416  0.0,result.get(),1);
417 
418  for (i = 0; i < nInteriorDofs; i++)
419  {
420  outarray[ mapArray[i] ] = result[i];
421  }
422  }
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:15
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:762
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:703
LibUtilities::NekManager< StdMatrixKey, DNekBlkMat, StdMatrixKey::opLess > m_stdStaticCondMatrixManager
virtual int v_NumBndryCoeffs() const
Definition: StdTriExp.cpp:793
virtual void v_GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdTriExp.cpp:998
std::shared_ptr< StdSegExp > StdSegExpSharedPtr
Definition: StdSegExp.h:46
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:372

References 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_GeneralMatrixOp_MatOp()

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

Reimplemented in Nektar::LocalRegions::TriExp.

Definition at line 1538 of file StdTriExp.cpp.

1542  {
1544 
1545  if(inarray.get() == outarray.get())
1546  {
1547  Array<OneD,NekDouble> tmp(m_ncoeffs);
1548  Vmath::Vcopy(m_ncoeffs,inarray.get(),1,tmp.get(),1);
1549 
1550  Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, 1.0, mat->GetPtr().get(),
1551  m_ncoeffs, tmp.get(), 1, 0.0, outarray.get(), 1);
1552  }
1553  else
1554  {
1555  Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, 1.0, mat->GetPtr().get(),
1556  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
1557  }
1558  }
LibUtilities::NekManager< StdMatrixKey, DNekMat, StdMatrixKey::opLess > m_stdMatrixManager

References Blas::Dgemv(), Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::StdRegions::StdExpansion::m_stdMatrixManager, and Vmath::Vcopy().

◆ v_GenMatrix()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 1246 of file StdTriExp.cpp.

1247  {
1248 
1249  MatrixType mtype = mkey.GetMatrixType();
1250 
1251  DNekMatSharedPtr Mat;
1252 
1253  switch(mtype)
1254  {
1256  {
1257  int nq0, nq1, nq;
1258 
1259  nq0 = m_base[0]->GetNumPoints();
1260  nq1 = m_base[1]->GetNumPoints();
1261 
1262  // take definition from key
1263  if(mkey.ConstFactorExists(eFactorConst))
1264  {
1265  nq = (int) mkey.GetConstFactor(eFactorConst);
1266  }
1267  else
1268  {
1269  nq = max(nq0,nq1);
1270  }
1271 
1272  int neq = LibUtilities::StdTriData::
1274  Array<OneD, Array<OneD, NekDouble> > coords(neq);
1275  Array<OneD, NekDouble> coll (2);
1276  Array<OneD, DNekMatSharedPtr> I (2);
1277  Array<OneD, NekDouble> tmp (nq0);
1278 
1279 
1280  Mat = MemoryManager<DNekMat>::AllocateSharedPtr(neq,nq0*nq1);
1281  int cnt = 0;
1282 
1283  for(int i = 0; i < nq; ++i)
1284  {
1285  for(int j = 0; j < nq-i; ++j,++cnt)
1286  {
1287  coords[cnt] = Array<OneD, NekDouble>(2);
1288  coords[cnt][0] = -1.0 + 2*j/(NekDouble)(nq-1);
1289  coords[cnt][1] = -1.0 + 2*i/(NekDouble)(nq-1);
1290  }
1291  }
1292 
1293  for(int i = 0; i < neq; ++i)
1294  {
1295  LocCoordToLocCollapsed(coords[i],coll);
1296 
1297  I[0] = m_base[0]->GetI(coll);
1298  I[1] = m_base[1]->GetI(coll+1);
1299 
1300  // interpolate first coordinate direction
1301  for (int j = 0; j < nq1; ++j)
1302  {
1303  NekDouble fac = (I[1]->GetPtr())[j];
1304  Vmath::Smul(nq0, fac, I[0]->GetPtr(), 1, tmp, 1);
1305 
1306  Vmath::Vcopy(nq0, &tmp[0], 1,
1307  Mat->GetRawPtr() + j*nq0*neq + i, neq);
1308  }
1309 
1310  }
1311  break;
1312  }
1313  default:
1314  {
1316  break;
1317  }
1318  }
1319 
1320  return Mat;
1321  }
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
Definition: StdExpansion.h:982
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:225

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp.

Definition at line 1028 of file StdTriExp.cpp.

1029  {
1032  "Expansion not of expected type");
1033  int i;
1034  int cnt;
1035  int nummodes0, nummodes1;
1036  int value;
1037 
1038  if (outarray.size()!=NumBndryCoeffs())
1039  {
1040  outarray = Array<OneD, unsigned int>(NumBndryCoeffs());
1041  }
1042 
1043  nummodes0 = m_base[0]->GetNumModes();
1044  nummodes1 = m_base[1]->GetNumModes();
1045 
1046  value = 2*nummodes1-1;
1047  for(i = 0; i < value; i++)
1048  {
1049  outarray[i]=i;
1050  }
1051  cnt = value;
1052 
1053  for(i = 0; i < nummodes0-2; i++)
1054  {
1055  outarray[cnt++]=value;
1056  value += nummodes1-2-i;
1057  }
1058  }

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 853 of file StdTriExp.cpp.

856  {
857  boost::ignore_unused(coords_2);
858 
859  Array<OneD, const NekDouble> z0 = m_base[0]->GetZ();
860  Array<OneD, const NekDouble> z1 = m_base[1]->GetZ();
861  int nq0 = GetNumPoints(0);
862  int nq1 = GetNumPoints(1);
863  int i,j;
864 
865  for(i = 0; i < nq1; ++i)
866  {
867  for(j = 0; j < nq0; ++j)
868  {
869  coords_0[i*nq0+j] = (1+z0[j])*(1-z1[i])/2.0 - 1.0;
870  }
871  Vmath::Fill(nq0,z1[i],&coords_1[0] + i*nq0,1);
872  }
873  }

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp.

Definition at line 998 of file StdTriExp.cpp.

999  {
1002  "Expansion not of a proper type");
1003 
1004  int i,j;
1005  int cnt = 0;
1006  int nummodes0, nummodes1;
1007  int startvalue;
1008  if(outarray.size()!=GetNcoeffs()-NumBndryCoeffs())
1009  {
1010  outarray = Array<OneD, unsigned int>(GetNcoeffs()-NumBndryCoeffs());
1011  }
1012 
1013  nummodes0 = m_base[0]->GetNumModes();
1014  nummodes1 = m_base[1]->GetNumModes();
1015 
1016  startvalue = 2*nummodes1;
1017 
1018  for(i = 0; i < nummodes0-2; i++)
1019  {
1020  for(j = 0; j < nummodes1-3-i; j++)
1021  {
1022  outarray[cnt++]=startvalue+j;
1023  }
1024  startvalue+=nummodes1-2-i;
1025  }
1026  }
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:124

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

◆ v_GetNtraces()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 783 of file StdTriExp.cpp.

784  {
785  return 3;
786  }

◆ v_GetNverts()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 778 of file StdTriExp.cpp.

779  {
780  return 3;
781  }

◆ v_GetSimplexEquiSpacedConnectivity()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1609 of file StdTriExp.cpp.

1612  {
1613  boost::ignore_unused(standard);
1614 
1615  int np1 = m_base[0]->GetNumPoints();
1616  int np2 = m_base[1]->GetNumPoints();
1617  int np = max(np1,np2);
1618 
1619  conn = Array<OneD, int>(3*(np-1)*(np-1));
1620 
1621  int row = 0;
1622  int rowp1 = 0;
1623  int cnt = 0;
1624  for(int i = 0; i < np-1; ++i)
1625  {
1626  rowp1 += np-i;
1627  for(int j = 0; j < np-i-2; ++j)
1628  {
1629  conn[cnt++] = row +j;
1630  conn[cnt++] = row +j+1;
1631  conn[cnt++] = rowp1 +j;
1632 
1633  conn[cnt++] = rowp1 +j+1;
1634  conn[cnt++] = rowp1 +j;
1635  conn[cnt++] = row +j+1;
1636  }
1637 
1638  conn[cnt++] = row +np-i-2;
1639  conn[cnt++] = row +np-i-1;
1640  conn[cnt++] = rowp1+np-i-2;
1641 
1642  row += np-i;
1643  }
1644  }

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

◆ v_GetTraceBasisKey()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 881 of file StdTriExp.cpp.

883  {
884  boost::ignore_unused(j);
885  ASSERTL2(i >= 0 && i <= 2, "edge id is out of range");
886 
887  if (i == 0)
888  {
889  return GetBasis(0)->GetBasisKey();
890  }
891  else
892  {
893  switch(m_base[1]->GetBasisType())
894  {
896  {
897  switch(m_base[1]->GetPointsType())
898  {
900  {
901  LibUtilities::PointsKey pkey(
902  m_base[1]->GetBasisKey().GetPointsKey().
903  GetNumPoints()+1,
905  return LibUtilities::BasisKey(
907  m_base[1]->GetNumModes(),pkey);
908  break;
909  }
910 
911  default:
913  "unexpected points distribution");
914  break;
915  }
916  break;
917  }
918  default:
920  "Information not available to set edge key");
921  break;
922  }
923  }
925  }
#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:111
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:208
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
@ eGaussRadauMAlpha1Beta0
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:58

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

◆ v_GetTraceInteriorToElementMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp.

Definition at line 1172 of file StdTriExp.cpp.

1177  {
1180  "Mapping not defined for this type of basis");
1181  int i;
1182  const int nummodes1 = m_base[1]->GetNumModes();
1183  const int nEdgeIntCoeffs = GetTraceNcoeffs(eid)-2;
1184 
1185  if(maparray.size() != nEdgeIntCoeffs)
1186  {
1187  maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
1188  }
1189 
1190  if(signarray.size() != nEdgeIntCoeffs)
1191  {
1192  signarray = Array<OneD, int>(nEdgeIntCoeffs,1);
1193  }
1194  else
1195  {
1196  fill( signarray.get() , signarray.get()+nEdgeIntCoeffs, 1 );
1197  }
1198 
1199  switch(eid)
1200  {
1201  case 0:
1202  {
1203  int cnt = 2*nummodes1 - 1;
1204  for(i = 0; i < nEdgeIntCoeffs; cnt+=nummodes1-2-i, ++i)
1205  {
1206  maparray[i] = cnt;
1207  }
1208  break;
1209  }
1210  case 1:
1211  {
1212  for(i = 0; i < nEdgeIntCoeffs; i++)
1213  {
1214  maparray[i] = nummodes1+1+i;
1215  }
1216  break;
1217  }
1218  case 2:
1219  {
1220  for(i = 0; i < nEdgeIntCoeffs; i++)
1221  {
1222  maparray[i] = 2+i;
1223  }
1224  break;
1225  }
1226  default:
1227  {
1228  ASSERTL0(false,"eid must be between 0 and 2");
1229  break;
1230  }
1231  }
1232 
1233  if(edgeOrient==eBackwards)
1234  {
1235  for(i = 1; i < nEdgeIntCoeffs; i+=2)
1236  {
1237  signarray[i] = -1;
1238  }
1239  }
1240  }
int GetTraceNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th trace.
Definition: StdExpansion.h:265

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 813 of file StdTriExp.cpp.

814  {
815  ASSERTL2(i >= 0 && i <= 2, "edge id is out of range");
816 
817  if (i == 0)
818  {
819  return GetBasisNumModes(0);
820  }
821  else
822  {
823  return GetBasisNumModes(1);
824  }
825  }
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:171

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

◆ v_GetTraceNumPoints()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 827 of file StdTriExp.cpp.

828  {
829  ASSERTL2((i >= 0)&&(i <= 2),"edge id is out of range");
830 
831  if (i == 0)
832  {
833  return GetNumPoints(0);
834  }
835  else
836  {
837  return GetNumPoints(1);
838  }
839  }

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

◆ v_GetTraceToElementMap()

void Nektar::StdRegions::StdTriExp::v_GetTraceToElementMap ( const int  eid,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray,
Orientation  edgeOrient = eForwards,
int  P = -1,
int  Q = -1 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp.

Definition at line 1060 of file StdTriExp.cpp.

1066  {
1067  boost::ignore_unused(Q);
1068 
1071  "Mapping not defined for this type of basis");
1072 
1073  int i;
1074  int numModes=0;
1075  int order0 = m_base[0]->GetNumModes();
1076  int order1 = m_base[1]->GetNumModes();
1077 
1078  switch (eid)
1079  {
1080  case 0:
1081  numModes = order0;
1082  break;
1083  case 1:
1084  case 2:
1085  numModes = order1;
1086  break;
1087  default:
1088  ASSERTL0(false,"eid must be between 0 and 2");
1089  }
1090 
1091  bool checkForZeroedModes = false;
1092  if (P == -1)
1093  {
1094  P = numModes;
1095  }
1096  else if(P != numModes)
1097  {
1098  checkForZeroedModes = true;
1099  }
1100 
1101 
1102  if(maparray.size() != P)
1103  {
1104  maparray = Array<OneD, unsigned int>(P);
1105  }
1106 
1107  if(signarray.size() != P)
1108  {
1109  signarray = Array<OneD, int>(P,1);
1110  }
1111  else
1112  {
1113  fill(signarray.get() , signarray.get()+P, 1);
1114  }
1115 
1116  switch(eid)
1117  {
1118  case 0:
1119  {
1120  int cnt = 0;
1121  for(i = 0; i < P; cnt+=order1-i, ++i)
1122  {
1123  maparray[i] = cnt;
1124  }
1125  break;
1126  }
1127  case 1:
1128  {
1129  maparray[0] = order1;
1130  maparray[1] = 1;
1131  for(i = 2; i < P; i++)
1132  {
1133  maparray[i] = order1-1+i;
1134  }
1135  break;
1136  }
1137  case 2:
1138  {
1139  for(i = 0; i < P; i++)
1140  {
1141  maparray[i] = i;
1142  }
1143  break;
1144  }
1145  default:
1146  ASSERTL0(false,"eid must be between 0 and 2");
1147  break;
1148  }
1149 
1150  if(edgeOrient==eBackwards)
1151  {
1152  swap( maparray[0] , maparray[1] );
1153 
1154  for(i = 3; i < P; i+=2)
1155  {
1156  signarray[i] = -1;
1157  }
1158  }
1159 
1160  if (checkForZeroedModes)
1161  {
1162  // Zero signmap and set maparray to zero if
1163  // elemental modes are not as large as face modes
1164  for (int j = numModes; j < P; j++)
1165  {
1166  signarray[j] = 0.0;
1167  maparray[j] = maparray[0];
1168  }
1169  }
1170  }
P
Definition: main.py:133

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

◆ v_GetVertexMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp.

Definition at line 934 of file StdTriExp.cpp.

935  {
936  ASSERTL1(
939  "Mapping not defined for this type of basis");
940 
941  int localDOF = 0;
942  if(useCoeffPacking == true)
943  {
944  switch(localVertexId)
945  {
946  case 0:
947  {
948  localDOF = 0;
949  break;
950  }
951  case 1:
952  {
953  localDOF = 1;
954  break;
955  }
956  case 2:
957  {
958  localDOF = m_base[1]->GetNumModes();
959  break;
960  }
961  default:
962  {
963  ASSERTL0(false,"eid must be between 0 and 2");
964  break;
965  }
966  }
967  }
968  else // follow book format for vertex indexing.
969  {
970  switch(localVertexId)
971  {
972  case 0:
973  {
974  localDOF = 0;
975  break;
976  }
977  case 1:
978  {
979  localDOF = m_base[1]->GetNumModes();
980  break;
981  }
982  case 2:
983  {
984  localDOF = 1;
985  break;
986  }
987  default:
988  {
989  ASSERTL0(false,"eid must be between 0 and 2");
990  break;
991  }
992  }
993  }
994 
995  return localDOF;
996  }

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 1369 of file StdTriExp.cpp.

1373  {
1374  StdTriExp::v_HelmholtzMatrixOp_MatFree(inarray,outarray,mkey);
1375  }
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)

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

◆ v_Integral()

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

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 84 of file StdTriExp.cpp.

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

References Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, 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 
)
protectedvirtual

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::StdRegions::StdNodalTriExp, Nektar::LocalRegions::TriExp, and Nektar::LocalRegions::NodalTriExp.

Definition at line 466 of file StdTriExp.cpp.

469  {
470  StdTriExp::v_IProductWRTBase_SumFac(inarray,outarray);
471  }
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
Definition: StdTriExp.cpp:485

References v_IProductWRTBase_SumFac().

Referenced by v_FwdTrans(), and v_FwdTrans_BndConstrained().

◆ v_IProductWRTBase_MatOp()

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

Reimplemented in Nektar::LocalRegions::TriExp.

Definition at line 473 of file StdTriExp.cpp.

476  {
477  int nq = GetTotPoints();
478  StdMatrixKey iprodmatkey(eIProductWRTBase,DetShapeType(),*this);
479  DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
480 
481  Blas::Dgemv('N',m_ncoeffs,nq,1.0,iprodmat->GetPtr().get(),
482  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
483  }
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:134

References Nektar::StdRegions::StdExpansion::DetShapeType(), Blas::Dgemv(), Nektar::StdRegions::eIProductWRTBase, Nektar::StdRegions::StdExpansion::GetStdMatrix(), Nektar::StdRegions::StdExpansion::GetTotPoints(), and Nektar::StdRegions::StdExpansion::m_ncoeffs.

◆ v_IProductWRTBase_SumFac()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 485 of file StdTriExp.cpp.

489  {
490  int nquad0 = m_base[0]->GetNumPoints();
491  int nquad1 = m_base[1]->GetNumPoints();
492  int order0 = m_base[0]->GetNumModes();
493 
494  if(multiplybyweights)
495  {
496  Array<OneD,NekDouble> tmp(nquad0*nquad1+nquad1*order0);
497  Array<OneD,NekDouble> wsp(tmp+nquad0*nquad1);
498 
499  // multiply by integration constants
500  MultiplyByQuadratureMetric(inarray,tmp);
501  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
502  m_base[1]->GetBdata(),
503  tmp,outarray,wsp);
504  }
505  else
506  {
507  Array<OneD,NekDouble> wsp(nquad1*order0);
508  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
509  m_base[1]->GetBdata(),
510  inarray,outarray,wsp);
511  }
512  }
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:733

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

Implements Nektar::StdRegions::StdExpansion2D.

Definition at line 514 of file StdTriExp.cpp.

522  {
523  boost::ignore_unused(doCheckCollDir0, doCheckCollDir1);
524 
525  int i;
526  int mode;
527  int nquad0 = m_base[0]->GetNumPoints();
528  int nquad1 = m_base[1]->GetNumPoints();
529  int nmodes0 = m_base[0]->GetNumModes();
530  int nmodes1 = m_base[1]->GetNumModes();
531 
532  ASSERTL1(wsp.size() >= nquad1*nmodes0,
533  "Workspace size is not sufficient");
534 
535  Blas::Dgemm('T','N',nquad1,nmodes0,nquad0,1.0,inarray.get(),nquad0,
536  base0.get(),nquad0,0.0,wsp.get(),nquad1);
537 
538  // Inner product with respect to 'b' direction
539  for (mode=i=0; i < nmodes0; ++i)
540  {
541  Blas::Dgemv('T',nquad1,nmodes1-i,1.0, base1.get()+mode*nquad1,
542  nquad1,wsp.get()+i*nquad1,1, 0.0,
543  outarray.get() + mode,1);
544  mode += nmodes1 - i;
545  }
546 
547  // fix for modified basis by splitting top vertex mode
549  {
550  outarray[1] += Blas::Ddot(nquad1,base1.get()+nquad1,1,
551  wsp.get()+nquad1,1);
552  }
553  }
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:197

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 555 of file StdTriExp.cpp.

559  {
560  StdTriExp::v_IProductWRTDerivBase_SumFac(dir,inarray,outarray);
561  }
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdTriExp.cpp:597

References v_IProductWRTDerivBase_SumFac().

◆ v_IProductWRTDerivBase_MatOp()

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

Reimplemented in Nektar::LocalRegions::TriExp.

Definition at line 563 of file StdTriExp.cpp.

567  {
568  int nq = GetTotPoints();
570 
571  switch(dir)
572  {
573  case 0:
574  {
575  mtype = eIProductWRTDerivBase0;
576  break;
577  }
578  case 1:
579  {
580  mtype = eIProductWRTDerivBase1;
581  break;
582  }
583  default:
584  {
585  ASSERTL1(false,"input dir is out of range");
586  break;
587  }
588  }
589 
590  StdMatrixKey iprodmatkey(mtype,DetShapeType(),*this);
591  DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
592 
593  Blas::Dgemv('N',m_ncoeffs,nq,1.0,iprodmat->GetPtr().get(),
594  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
595  }

References ASSERTL1, Nektar::StdRegions::StdExpansion::DetShapeType(), Blas::Dgemv(), Nektar::StdRegions::eIProductWRTDerivBase0, Nektar::StdRegions::eIProductWRTDerivBase1, Nektar::StdRegions::StdExpansion::GetStdMatrix(), Nektar::StdRegions::StdExpansion::GetTotPoints(), and Nektar::StdRegions::StdExpansion::m_ncoeffs.

◆ v_IProductWRTDerivBase_SumFac()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 597 of file StdTriExp.cpp.

601  {
602  int i;
603  int nquad0 = m_base[0]->GetNumPoints();
604  int nquad1 = m_base[1]->GetNumPoints();
605  int nqtot = nquad0*nquad1;
606  int nmodes0 = m_base[0]->GetNumModes();
607  int wspsize = max(max(nqtot,m_ncoeffs),nquad1*nmodes0);
608 
609  Array<OneD, NekDouble> gfac0(2*wspsize);
610  Array<OneD, NekDouble> tmp0 (gfac0+wspsize);
611 
612  const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
613 
614  // set up geometric factor: 2/(1-z1)
615  for(i = 0; i < nquad1; ++i)
616  {
617  gfac0[i] = 2.0/(1-z1[i]);
618  }
619 
620  for(i = 0; i < nquad1; ++i)
621  {
622  Vmath::Smul(nquad0,gfac0[i],&inarray[0]+i*nquad0,1,
623  &tmp0[0]+i*nquad0,1);
624  }
625 
626  MultiplyByQuadratureMetric(tmp0,tmp0);
627 
628  switch(dir)
629  {
630  case 0:
631  {
632  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
633  m_base[1]->GetBdata(),
634  tmp0,outarray,gfac0);
635  break;
636  }
637  case 1:
638  {
639  Array<OneD, NekDouble> tmp3(m_ncoeffs);
640  const Array<OneD, const NekDouble>& z0 = m_base[0]->GetZ();
641 
642  for (i = 0; i < nquad0; ++i)
643  {
644  gfac0[i] = 0.5*(1+z0[i]);
645  }
646 
647  for (i = 0; i < nquad1; ++i)
648  {
649  Vmath::Vmul(nquad0,&gfac0[0],1,&tmp0[0]+i*nquad0,1,
650  &tmp0[0]+i*nquad0,1);
651  }
652 
653  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
654  m_base[1]->GetBdata(),
655  tmp0,tmp3,gfac0);
656 
657  MultiplyByQuadratureMetric(inarray,tmp0);
658  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
659  m_base[1]->GetDbdata(),
660  tmp0,outarray,gfac0);
661  Vmath::Vadd(m_ncoeffs,&tmp3[0],1,&outarray[0],1,
662  &outarray[0],1);
663  break;
664  }
665  default:
666  {
667  ASSERTL1(false, "input dir is out of range");
668  break;
669  }
670  }
671  }
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:322

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 875 of file StdTriExp.cpp.

876  {
877  return m_base[0]->GetBasisType() == LibUtilities::eModified_A &&
878  m_base[1]->GetBasisType() == LibUtilities::eModified_B;
879  }

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 1341 of file StdTriExp.cpp.

1345  {
1346  StdTriExp::v_LaplacianMatrixOp_MatFree(inarray,outarray,mkey);
1347  }
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 1349 of file StdTriExp.cpp.

1355  {
1357  k1,k2,inarray,outarray,mkey);
1358  }
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 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 696 of file StdTriExp.cpp.

698  {
699  xi[0] = (1.0 + eta[0]) * (1.0 - eta[1]) * 0.5 - 1.0;
700  xi[1] = eta[1];
701  }

◆ v_LocCoordToLocCollapsed()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 677 of file StdTriExp.cpp.

679  {
680  NekDouble d1 = 1.-xi[1];
681  if(fabs(d1) < NekConstants::kNekZeroTol)
682  {
683  if(d1>=0.)
684  {
686  }
687  else
688  {
690  }
691  }
692  eta[0] = 2.*(1.+xi[0])/d1-1.0;
693  eta[1] = xi[1];
694  }
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 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 1333 of file StdTriExp.cpp.

1337  {
1338  StdExpansion::MassMatrixOp_MatFree(inarray,outarray,mkey);
1339  }
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 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1564 of file StdTriExp.cpp.

1567  {
1568  int i;
1569  int nquad0 = m_base[0]->GetNumPoints();
1570  int nquad1 = m_base[1]->GetNumPoints();
1571 
1572  const Array<OneD, const NekDouble>& w0 = m_base[0]->GetW();
1573  const Array<OneD, const NekDouble>& w1 = m_base[1]->GetW();
1574  const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
1575 
1576  // multiply by integration constants
1577  for(i = 0; i < nquad1; ++i)
1578  {
1579  Vmath::Vmul(nquad0,inarray.get()+i*nquad0,1,
1580  w0.get(),1, outarray.get()+i*nquad0,1);
1581  }
1582 
1583  switch(m_base[1]->GetPointsType())
1584  {
1585  // Legendre inner product
1588  for(i = 0; i < nquad1; ++i)
1589  {
1590  Blas::Dscal(nquad0,0.5*(1-z1[i])*w1[i],
1591  outarray.get()+i*nquad0,1);
1592  }
1593  break;
1594 
1595  // (1,0) Jacobi Inner product
1597  for(i = 0; i < nquad1; ++i)
1598  {
1599  Blas::Dscal(nquad0,0.5*w1[i],outarray.get()+i*nquad0,1);
1600  }
1601  break;
1602 
1603  default:
1604  ASSERTL0(false, "Unsupported quadrature points type.");
1605  break;
1606  }
1607  }
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition: Blas.hpp:182
@ ePolyEvenlySpaced
1D Evenly-spaced points using Lagrange polynomial
Definition: PointsType.h:64

References ASSERTL0, Blas::Dscal(), Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, Nektar::LibUtilities::ePolyEvenlySpaced, Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::StdRegions::StdExpansion::m_base, and Vmath::Vmul().

◆ v_NumBndryCoeffs()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp.

Definition at line 793 of file StdTriExp.cpp.

794  {
796  "BasisType is not a boundary interior form");
798  "BasisType is not a boundary interior form");
799 
800  return 3 + (GetBasisNumModes(0)-2) + 2*(GetBasisNumModes(1)-2);
801  }

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

Referenced by v_FwdTrans_BndConstrained().

◆ v_NumDGBndryCoeffs()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 803 of file StdTriExp.cpp.

804  {
806  "BasisType is not a boundary interior form");
808  "BasisType is not a boundary interior form");
809 
810  return GetBasisNumModes(0) + 2*GetBasisNumModes(1);
811  }

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

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 132 of file StdTriExp.cpp.

137  {
138  boost::ignore_unused(out_d2);
139 
140  int i;
141  int nquad0 = m_base[0]->GetNumPoints();
142  int nquad1 = m_base[1]->GetNumPoints();
143  Array<OneD, NekDouble> wsp(std::max(nquad0, nquad1));
144 
145  const Array<OneD, const NekDouble>& z0 = m_base[0]->GetZ();
146  const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
147 
148  // set up geometric factor: 2/(1-z1)
149  Vmath::Sadd(nquad1, -1.0, z1, 1, wsp, 1);
150  Vmath::Sdiv(nquad1, -2.0, wsp, 1, wsp, 1);
151 
152  if (out_d0.size() > 0)
153  {
154  PhysTensorDeriv(inarray, out_d0, out_d1);
155 
156  for (i = 0; i < nquad1; ++i)
157  {
158  Blas::Dscal(nquad0,wsp[i],&out_d0[0]+i*nquad0,1);
159  }
160 
161  // if no d1 required do not need to calculate both deriv
162  if (out_d1.size() > 0)
163  {
164  // set up geometric factor: (1_z0)/(1-z1)
165  Vmath::Sadd(nquad0, 1.0, z0, 1, wsp, 1);
166  Vmath::Smul(nquad0, 0.5, wsp, 1, wsp, 1);
167 
168  for (i = 0; i < nquad1; ++i)
169  {
170  Vmath::Vvtvp(nquad0,&wsp[0],1,&out_d0[0]+i*nquad0,
171  1,&out_d1[0]+i*nquad0,
172  1,&out_d1[0]+i*nquad0,1);
173  }
174  }
175  }
176  else if (out_d1.size() > 0)
177  {
178  Array<OneD, NekDouble> diff0(nquad0*nquad1);
179  PhysTensorDeriv(inarray, diff0, out_d1);
180 
181  for (i = 0; i < nquad1; ++i)
182  {
183  Blas::Dscal(nquad0,wsp[i],&diff0[0]+i*nquad0,1);
184  }
185 
186  Vmath::Sadd(nquad0, 1.0, z0, 1, wsp, 1);
187  Vmath::Smul(nquad0, 0.5, wsp, 1, wsp, 1);
188 
189  for (i = 0; i < nquad1; ++i)
190  {
191  Vmath::Vvtvp(nquad0,&wsp[0],1,&diff0[0]+i*nquad0,
192  1,&out_d1[0]+i*nquad0,
193  1,&out_d1[0]+i*nquad0,1);
194  }
195  }
196  }
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:513
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:291
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha - x.
Definition: Vmath.cpp:341

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

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

See also
StdRegions::StdExpansion::PhysDeriv

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 198 of file StdTriExp.cpp.

202  {
203  switch(dir)
204  {
205  case 0:
206  {
207  v_PhysDeriv(inarray, outarray, NullNekDouble1DArray);
208  break;
209  }
210  case 1:
211  {
212  v_PhysDeriv(inarray, NullNekDouble1DArray, outarray);
213  break;
214  }
215  default:
216  {
217  ASSERTL1(false,"input dir is out of range");
218  break;
219  }
220  }
221  }
virtual 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)
Calculate the derivative of the physical points.
Definition: StdTriExp.cpp:132
static Array< OneD, NekDouble > NullNekDouble1DArray

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

◆ v_PhysEvaluateBasis()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 751 of file StdTriExp.cpp.

754  {
755  Array<OneD, NekDouble> coll(2);
756  LocCoordToLocCollapsed(coords, coll);
757 
758  // From mode we need to determine mode0 and mode1 in the (p,q)
759  // direction. mode1 can be directly inferred from mode.
760  const int nm1 = m_base[1]->GetNumModes();
761  const double c = 1 + 2*nm1;
762  const int mode0 = floor(0.5*(c - sqrt(c*c - 8*mode)));
763 
764  if (mode == 1 &&
766  {
767  // Account for collapsed vertex.
768  return StdExpansion::BaryEvaluateBasis<1>(coll[1], 1);
769  }
770  else
771  {
772  return
773  StdExpansion::BaryEvaluateBasis<0>(coll[0], mode0) *
774  StdExpansion::BaryEvaluateBasis<1>(coll[1], mode);
775  }
776  }
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:267

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp.

Definition at line 1481 of file StdTriExp.cpp.

1485  {
1486  int n_coeffs = inarray.size();
1487  int nquad0 = m_base[0]->GetNumPoints();
1488  int nquad1 = m_base[1]->GetNumPoints();
1489  Array<OneD, NekDouble> coeff(n_coeffs);
1490  Array<OneD, NekDouble> coeff_tmp(n_coeffs,0.0);
1491  Array<OneD, NekDouble> tmp;
1492  Array<OneD, NekDouble> tmp2;
1493  int nqtot = nquad0*nquad1;
1494  Array<OneD, NekDouble> phys_tmp(nqtot,0.0);
1495 
1496  int nmodes0 = m_base[0]->GetNumModes();
1497  int nmodes1 = m_base[1]->GetNumModes();
1498  int numMin2 = nmodes0;
1499  int i;
1500 
1501  const LibUtilities::PointsKey Pkey0(
1503  const LibUtilities::PointsKey Pkey1(
1505 
1506  LibUtilities::BasisKey b0(m_base[0]->GetBasisType(),nmodes0,Pkey0);
1507  LibUtilities::BasisKey b1(m_base[1]->GetBasisType(),nmodes1,Pkey1);
1508 
1509  LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A,nmodes0,Pkey0);
1510  LibUtilities::BasisKey bortho1(LibUtilities::eOrtho_B,nmodes1,Pkey1);
1511 
1512  StdRegions::StdTriExpSharedPtr m_OrthoTriExp;
1514 
1516  ::AllocateSharedPtr(b0, b1);
1517  m_OrthoTriExp = MemoryManager<StdRegions::StdTriExp>
1518  ::AllocateSharedPtr(bortho0, bortho1);
1519 
1520  m_TriExp ->BwdTrans(inarray,phys_tmp);
1521  m_OrthoTriExp->FwdTrans(phys_tmp, coeff);
1522 
1523  for (i = 0; i < n_coeffs; i++)
1524  {
1525  if (i == numMin)
1526  {
1527  coeff[i] = 0.0;
1528  numMin += numMin2 - 1;
1529  numMin2 -= 1.0;
1530  }
1531  }
1532 
1533  m_OrthoTriExp->BwdTrans(coeff,phys_tmp);
1534  m_TriExp ->FwdTrans(phys_tmp, outarray);
1535 
1536  }
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:45
std::shared_ptr< StdTriExp > StdTriExpSharedPtr
Definition: StdTriExp.h:272

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::NodalTriExp.

Definition at line 223 of file StdTriExp.cpp.

228  {
229  boost::ignore_unused(out_d2);
230  StdTriExp::v_PhysDeriv(inarray, out_d0, out_d1);
231  }

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 233 of file StdTriExp.cpp.

237  {
238  StdTriExp::v_PhysDeriv(dir,inarray,outarray);
239  }

References v_PhysDeriv().

◆ v_SVVLaplacianFilter()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp.

Definition at line 1378 of file StdTriExp.cpp.

1380  {
1381  int qa = m_base[0]->GetNumPoints();
1382  int qb = m_base[1]->GetNumPoints();
1383  int nmodes_a = m_base[0]->GetNumModes();
1384  int nmodes_b = m_base[1]->GetNumModes();
1385 
1386  // Declare orthogonal basis.
1387  LibUtilities::PointsKey pa(qa,m_base[0]->GetPointsType());
1388  LibUtilities::PointsKey pb(qb,m_base[1]->GetPointsType());
1389 
1390  LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A,nmodes_a,pa);
1391  LibUtilities::BasisKey Bb(LibUtilities::eOrtho_B,nmodes_b,pb);
1392  StdTriExp OrthoExp(Ba,Bb);
1393 
1394  Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
1395 
1396  // project onto physical space.
1397  OrthoExp.FwdTrans(array,orthocoeffs);
1398 
1399  if(mkey.ConstFactorExists(eFactorSVVPowerKerDiffCoeff)) // Rodrigo's power kern
1400  {
1401  NekDouble cutoff = mkey.GetConstFactor(eFactorSVVCutoffRatio);
1402  NekDouble SvvDiffCoeff =
1403  mkey.GetConstFactor(eFactorSVVPowerKerDiffCoeff)*
1404  mkey.GetConstFactor(eFactorSVVDiffCoeff);
1405 
1406  int cnt = 0;
1407  for(int j = 0; j < nmodes_a; ++j)
1408  {
1409  for(int k = 0; k < nmodes_b-j; ++k, ++cnt)
1410  {
1411  NekDouble fac = std::max(
1412  pow((1.0*j)/(nmodes_a-1),cutoff*nmodes_a),
1413  pow((1.0*k)/(nmodes_b-1),cutoff*nmodes_b));
1414 
1415  orthocoeffs[cnt] *= (SvvDiffCoeff *fac);
1416  }
1417  }
1418  }
1419  else if(mkey.ConstFactorExists(eFactorSVVDGKerDiffCoeff)) // Rodrigo/mansoor's DG kernel
1420  {
1421  NekDouble SvvDiffCoeff =
1422  mkey.GetConstFactor(eFactorSVVDGKerDiffCoeff)*
1423  mkey.GetConstFactor(eFactorSVVDiffCoeff);
1424  int max_ab = max(nmodes_a-kSVVDGFiltermodesmin,
1425  nmodes_b-kSVVDGFiltermodesmin);
1426  max_ab = max(max_ab,0);
1427  max_ab = min(max_ab,kSVVDGFiltermodesmax-kSVVDGFiltermodesmin);
1428 
1429  int cnt = 0;
1430  for(int j = 0; j < nmodes_a; ++j)
1431  {
1432  for(int k = 0; k < nmodes_b-j; ++k, ++cnt)
1433  {
1434  int maxjk = max(j,k);
1435  maxjk = min(maxjk,kSVVDGFiltermodesmax-1);
1436 
1437  orthocoeffs[cnt] *= SvvDiffCoeff *
1438  kSVVDGFilter[max_ab][maxjk];
1439  }
1440  }
1441  }
1442  else
1443  {
1444  NekDouble SvvDiffCoeff =
1445  mkey.GetConstFactor(eFactorSVVDiffCoeff);
1446 
1447  int cutoff = (int) (mkey.GetConstFactor(eFactorSVVCutoffRatio)*
1448  min(nmodes_a,nmodes_b));
1449 
1450  NekDouble epsilon = 1.0;
1451  int nmodes = min(nmodes_a,nmodes_b);
1452 
1453  int cnt = 0;
1454 
1455  // apply SVV filter (JEL)
1456  for(int j = 0; j < nmodes_a; ++j)
1457  {
1458  for(int k = 0; k < nmodes_b-j; ++k)
1459  {
1460  if(j + k >= cutoff)
1461  {
1462  orthocoeffs[cnt] *= (SvvDiffCoeff
1463  *exp(-(j+k-nmodes)*(j+k-nmodes)
1464  /((NekDouble)((j+k-cutoff+epsilon)
1465  *(j+k-cutoff+epsilon)))));
1466  }
1467  else
1468  {
1469  orthocoeffs[cnt] *= 0.0;
1470  }
1471  cnt++;
1472  }
1473  }
1474 
1475  }
1476 
1477  // backward transform to physical space
1478  OrthoExp.BwdTrans(orthocoeffs,array);
1479  }
const int kSVVDGFiltermodesmin
Definition: StdRegions.hpp:388
const int kSVVDGFiltermodesmax
Definition: StdRegions.hpp:389
const NekDouble kSVVDGFilter[9][11]
Definition: StdRegions.hpp:391

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 1360 of file StdTriExp.cpp.

1365  {
1366  StdExpansion::WeakDerivMatrixOp_MatFree(i,inarray,outarray,mkey);
1367  }
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().