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 (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb)
 
 StdTriExp ()=default
 
 StdTriExp (const StdTriExp &T)=default
 
 ~StdTriExp () override=default
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion2D
 StdExpansion2D (int numcoeffs, const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb)
 
 StdExpansion2D ()=default
 
 StdExpansion2D (const StdExpansion2D &T)=default
 
 ~StdExpansion2D () override=default
 
void PhysTensorDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d0, Array< OneD, NekDouble > &outarray_d1)
 Calculate the 2D derivative in the local tensor/collapsed coordinate at the physical points. More...
 
NekDouble Integral (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &w0, const Array< OneD, const NekDouble > &w1)
 
NekDouble BaryTensorDeriv (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
 
void BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
 
void IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion
 StdExpansion ()
 Default Constructor. More...
 
 StdExpansion (const int numcoeffs, const int numbases, const LibUtilities::BasisKey &Ba=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bb=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bc=LibUtilities::NullBasisKey)
 Constructor. More...
 
 StdExpansion (const StdExpansion &T)
 Copy Constructor. More...
 
virtual ~StdExpansion ()
 Destructor. More...
 
int GetNumBases () const
 This function returns the number of 1D bases used in the expansion. More...
 
const Array< OneD, const LibUtilities::BasisSharedPtr > & GetBase () const
 This function gets the shared point to basis. More...
 
const LibUtilities::BasisSharedPtrGetBasis (int dir) const
 This function gets the shared point to basis in the dir direction. More...
 
int GetNcoeffs (void) const
 This function returns the total number of coefficients used in the expansion. More...
 
int GetTotPoints () const
 This function returns the total number of quadrature points used in the element. More...
 
LibUtilities::BasisType GetBasisType (const int dir) const
 This function returns the type of basis used in the dir direction. More...
 
int GetBasisNumModes (const int dir) const
 This function returns the number of expansion modes in the dir direction. More...
 
int EvalBasisNumModesMax (void) const
 This function returns the maximum number of expansion modes over all local directions. More...
 
LibUtilities::PointsType GetPointsType (const int dir) const
 This function returns the type of quadrature points used in the dir direction. More...
 
int GetNumPoints (const int dir) const
 This function returns the number of quadrature points in the dir direction. More...
 
const Array< OneD, const NekDouble > & GetPoints (const int dir) const
 This function returns a pointer to the array containing the quadrature points in dir direction. More...
 
int GetNverts () const
 This function returns the number of vertices of the expansion domain. More...
 
int GetTraceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th trace. More...
 
int GetTraceIntNcoeffs (const int i) const
 
int GetTraceNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th trace. More...
 
const LibUtilities::BasisKey GetTraceBasisKey (const int i, int k=-1) const
 This function returns the basis key belonging to the i-th trace. More...
 
LibUtilities::PointsKey GetTracePointsKey (const int i, int k=-1) const
 This function returns the basis key belonging to the i-th trace. More...
 
int NumBndryCoeffs (void) const
 
int NumDGBndryCoeffs (void) const
 
const LibUtilities::PointsKey GetNodalPointsKey () const
 This function returns the type of expansion Nodal point type if defined. More...
 
int GetNtraces () const
 Returns the number of trace elements connected to this element. More...
 
LibUtilities::ShapeType DetShapeType () const
 This function returns the shape of the expansion domain. More...
 
std::shared_ptr< StdExpansionGetStdExp () const
 
std::shared_ptr< StdExpansionGetLinStdExp (void) const
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion () const
 
bool IsNodalNonTensorialExp ()
 
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space. More...
 
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space. More...
 
void FwdTransBndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble Integral (const Array< OneD, const NekDouble > &inarray)
 This function integrates the specified function over the domain. More...
 
void FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 This function fills the array outarray with the mode-th mode of the expansion. More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 this function calculates the inner product of a given function f with the different modes of the expansion More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check)
 
void IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
int GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void SetElmtId (const int id)
 Set the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2=NullNekDouble1DArray, Array< OneD, NekDouble > &coords_3=NullNekDouble1DArray)
 this function returns the physical coordinates of the quadrature points of the expansion More...
 
void GetCoord (const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
 given the coordinates of a point of the element in the local collapsed coordinate system, this function calculates the physical coordinates of the point More...
 
DNekMatSharedPtr GetStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr GetStdStaticCondMatrix (const StdMatrixKey &mkey)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
void DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
int CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
NekDouble StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
int GetCoordim ()
 
void GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
void GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
int GetVertexMap (const int localVertexId, bool useCoeffPacking=false)
 
void GetTraceToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
 
void GetTraceCoeffMap (const unsigned int traceid, Array< OneD, unsigned int > &maparray)
 
void GetElmtTraceToTraceMap (const unsigned int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
 
void GetTraceInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eForwards)
 
void GetTraceNumModes (const int tid, int &numModes0, int &numModes1, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2)
 
void MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
DNekMatSharedPtr CreateGeneralMatrix (const StdMatrixKey &mkey)
 this function generates the mass matrix \(\mathbf{M}[i][j] = \int \phi_i(\mathbf{x}) \phi_j(\mathbf{x}) d\mathbf{x}\) More...
 
void GeneralMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
void ExponentialFilter (Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
 
void LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
DNekMatSharedPtr GenMatrix (const StdMatrixKey &mkey)
 
void PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void PhysDeriv_s (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds)
 
void PhysDeriv_n (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn)
 
void PhysDirectionalDeriv (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &outarray)
 
void StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
 This function evaluates the first derivative of the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs, std::array< NekDouble, 6 > &secondOrderDerivs)
 
NekDouble PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluateBasis (const Array< OneD, const NekDouble > &coords, int mode)
 This function evaluates the basis function mode mode at a point coords of the domain. More...
 
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta. More...
 
void LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi)
 Convert local collapsed coordinates eta into local cartesian coordinate xi. More...
 
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
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

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

Additional Inherited Members

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

Detailed Description

Definition at line 43 of file StdTriExp.h.

Constructor & Destructor Documentation

◆ StdTriExp() [1/3]

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

Definition at line 41 of file StdTriExp.cpp.

44 Ba.GetNumModes(), Bb.GetNumModes()),
45 2, Ba, Bb),
47 Ba.GetNumModes(), Bb.GetNumModes()),
48 Ba, Bb)
49{
50 ASSERTL0(Ba.GetNumModes() <= Bb.GetNumModes(),
51 "order in 'a' direction is higher than order "
52 "in 'b' direction");
53}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
StdExpansion()
Default Constructor.
int getNumberOfCoefficients(int Na, int Nb)
Definition: ShapeType.hpp:109

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

◆ StdTriExp() [2/3]

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

◆ StdTriExp() [3/3]

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

◆ ~StdTriExp()

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

Member Function Documentation

◆ v_BwdTrans()

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

Backward tranform for triangular elements.

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 217 of file StdTriExp.cpp.

219{
220 v_BwdTrans_SumFac(inarray, outarray);
221}
void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdTriExp.cpp:223

References v_BwdTrans_SumFac().

◆ v_BwdTrans_SumFac()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::NodalTriExp.

Definition at line 223 of file StdTriExp.cpp.

225{
226 Array<OneD, NekDouble> wsp(m_base[1]->GetNumPoints() *
227 m_base[0]->GetNumModes());
228
229 BwdTrans_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetBdata(), inarray,
230 outarray, wsp);
231}
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:218
Array< OneD, LibUtilities::BasisSharedPtr > m_base

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

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

◆ v_BwdTrans_SumFacKernel()

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

Implements Nektar::StdRegions::StdExpansion2D.

Definition at line 233 of file StdTriExp.cpp.

240{
241 int i;
242 int mode;
243 int nquad0 = m_base[0]->GetNumPoints();
244 int nquad1 = m_base[1]->GetNumPoints();
245 int nmodes0 = m_base[0]->GetNumModes();
246 int nmodes1 = m_base[1]->GetNumModes();
247
248 ASSERTL1(wsp.size() >= nquad1 * nmodes0,
249 "Workspace size is not sufficient");
252 "Basis[1] is not of general tensor type");
253
254 for (i = mode = 0; i < nmodes0; ++i)
255 {
256 Blas::Dgemv('N', nquad1, nmodes1 - i, 1.0, base1.get() + mode * nquad1,
257 nquad1, &inarray[0] + mode, 1, 0.0, &wsp[0] + i * nquad1,
258 1);
259 mode += nmodes1 - i;
260 }
261
262 // fix for modified basis by splitting top vertex mode
264 {
265 Blas::Daxpy(nquad1, inarray[1], base1.get() + nquad1, 1,
266 &wsp[0] + nquad1, 1);
267 }
268
269 Blas::Dgemm('N', 'T', nquad0, nquad1, nmodes0, 1.0, base0.get(), nquad0,
270 &wsp[0], nquad1, 0.0, &outarray[0], nquad0);
271}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:265
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:156
static void Dgemv(const char &trans, const int &m, const int &n, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = alpha A x plus beta y where A[m x n].
Definition: Blas.hpp:211
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:383
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
Definition: Blas.hpp:135
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:49
@ eOrtho_B
Principle Orthogonal Functions .
Definition: BasisType.h:44
@ 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 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 794 of file StdTriExp.cpp.

796{
798 nummodes[modes_offset], nummodes[modes_offset + 1]);
799 modes_offset += 2;
800
801 return nmodes;
802}

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

◆ v_CreateStdMatrix()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 1295 of file StdTriExp.cpp.

1296{
1297 return v_GenMatrix(mkey);
1298}
DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override
Definition: StdTriExp.cpp:1221

References v_GenMatrix().

◆ v_DetShapeType()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 727 of file StdTriExp.cpp.

728{
730}

References Nektar::LibUtilities::eTriangle.

◆ v_FillMode()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 603 of file StdTriExp.cpp.

604{
605 int i, m;
606 int nquad0 = m_base[0]->GetNumPoints();
607 int nquad1 = m_base[1]->GetNumPoints();
608 int order0 = m_base[0]->GetNumModes();
609 int order1 = m_base[1]->GetNumModes();
610 int mode0 = 0;
611 Array<OneD, const NekDouble> base0 = m_base[0]->GetBdata();
612 Array<OneD, const NekDouble> base1 = m_base[1]->GetBdata();
613
614 ASSERTL2(mode <= m_ncoeffs, "calling argument mode is larger than "
615 "total expansion order");
616
617 m = order1;
618 for (i = 0; i < order0; ++i, m += order1 - i)
619 {
620 if (m > mode)
621 {
622 mode0 = i;
623 break;
624 }
625 }
626
627 // deal with top vertex mode in modified basis
628 if (mode == 1 && m_base[0]->GetBasisType() == LibUtilities::eModified_A)
629 {
630 Vmath::Fill(nquad0 * nquad1, 1.0, outarray, 1);
631 }
632 else
633 {
634 for (i = 0; i < nquad1; ++i)
635 {
636 Vmath::Vcopy(nquad0, (NekDouble *)(base0.get() + mode0 * nquad0), 1,
637 &outarray[0] + i * nquad0, 1);
638 }
639 }
640
641 for (i = 0; i < nquad0; ++i)
642 {
643 Vmath::Vmul(nquad1, (NekDouble *)(base1.get() + mode * nquad1), 1,
644 &outarray[0] + i, nquad0, &outarray[0] + i, nquad0);
645 }
646}
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.hpp:72
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.hpp:54
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

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

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

◆ v_FwdTrans()

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

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

See also
StdExpansion::FwdTrans

Implements Nektar::StdRegions::StdExpansion.

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

Definition at line 273 of file StdTriExp.cpp.

275{
276 v_IProductWRTBase(inarray, outarray);
277
278 // get Mass matrix inverse
279 StdMatrixKey masskey(eInvMass, DetShapeType(), *this);
280 DNekMatSharedPtr matsys = GetStdMatrix(masskey);
281
282 // copy inarray in case inarray == outarray
283 NekVector<NekDouble> in(m_ncoeffs, outarray, eCopy);
284 NekVector<NekDouble> out(m_ncoeffs, outarray, eWrapper);
285
286 out = (*matsys) * in;
287}
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:603
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:367
void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Calculate the inner product of inarray with respect to the basis B=base0[p]*base1[pq] and put into ou...
Definition: StdTriExp.cpp:420
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75

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

◆ v_FwdTransBndConstrained()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp.

Definition at line 289 of file StdTriExp.cpp.

292{
293 int i, j;
294 int npoints[2] = {m_base[0]->GetNumPoints(), m_base[1]->GetNumPoints()};
295 int nmodes[2] = {m_base[0]->GetNumModes(), m_base[1]->GetNumModes()};
296
297 fill(outarray.get(), outarray.get() + m_ncoeffs, 0.0);
298
299 Array<OneD, NekDouble> physEdge[3];
300 Array<OneD, NekDouble> coeffEdge[3];
301 for (i = 0; i < 3; i++)
302 {
303 physEdge[i] = Array<OneD, NekDouble>(npoints[i != 0]);
304 coeffEdge[i] = Array<OneD, NekDouble>(nmodes[i != 0]);
305 }
306
307 for (i = 0; i < npoints[0]; i++)
308 {
309 physEdge[0][i] = inarray[i];
310 }
311
312 for (i = 0; i < npoints[1]; i++)
313 {
314 physEdge[1][i] = inarray[npoints[0] - 1 + i * npoints[0]];
315 physEdge[2][i] =
316 inarray[(npoints[1] - 1) * npoints[0] - i * npoints[0]];
317 }
318
319 StdSegExpSharedPtr segexp[2] = {
321 m_base[0]->GetBasisKey()),
323 m_base[1]->GetBasisKey())};
324
325 Array<OneD, unsigned int> mapArray;
326 Array<OneD, int> signArray;
328
329 for (i = 0; i < 3; i++)
330 {
331 segexp[i != 0]->FwdTransBndConstrained(physEdge[i], coeffEdge[i]);
332
333 GetTraceToElementMap(i, mapArray, signArray);
334 for (j = 0; j < nmodes[i != 0]; j++)
335 {
336 sign = (NekDouble)signArray[j];
337 outarray[mapArray[j]] = sign * coeffEdge[i][j];
338 }
339 }
340
341 Array<OneD, NekDouble> tmp0(m_ncoeffs);
342 Array<OneD, NekDouble> tmp1(m_ncoeffs);
343
344 StdMatrixKey masskey(eMass, DetShapeType(), *this);
345 MassMatrixOp(outarray, tmp0, masskey);
346 v_IProductWRTBase(inarray, tmp1);
347
348 Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
349
350 // get Mass matrix inverse (only of interior DOF)
351 // use block (1,1) of the static condensed system
352 // note: this block alreay contains the inverse matrix
353 DNekMatSharedPtr matsys =
354 (m_stdStaticCondMatrixManager[masskey])->GetBlock(1, 1);
355
356 int nBoundaryDofs = v_NumBndryCoeffs();
357 int nInteriorDofs = m_ncoeffs - nBoundaryDofs;
358
359 Array<OneD, NekDouble> rhs(nInteriorDofs);
360 Array<OneD, NekDouble> result(nInteriorDofs);
361
362 v_GetInteriorMap(mapArray);
363
364 for (i = 0; i < nInteriorDofs; i++)
365 {
366 rhs[i] = tmp1[mapArray[i]];
367 }
368
369 Blas::Dgemv('N', nInteriorDofs, nInteriorDofs, 1.0, &(matsys->GetPtr())[0],
370 nInteriorDofs, rhs.get(), 1, 0.0, result.get(), 1);
371
372 for (i = 0; i < nInteriorDofs; i++)
373 {
374 outarray[mapArray[i]] = result[i];
375 }
376}
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:47
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:752
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:684
LibUtilities::NekManager< StdMatrixKey, DNekBlkMat, StdMatrixKey::opLess > m_stdStaticCondMatrixManager
void v_GetInteriorMap(Array< OneD, unsigned int > &outarray) override
Definition: StdTriExp.cpp:1032
int v_NumBndryCoeffs() const override
Definition: StdTriExp.cpp:732
std::shared_ptr< StdSegExp > StdSegExpSharedPtr
Definition: StdSegExp.h:222
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.hpp:220

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::StdRegions::StdExpansion::DetShapeType(), Blas::Dgemv(), Nektar::StdRegions::eMass, Nektar::StdRegions::StdExpansion::GetTraceToElementMap(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::StdRegions::StdExpansion::m_stdStaticCondMatrixManager, Nektar::StdRegions::StdExpansion::MassMatrixOp(), sign, v_GetInteriorMap(), v_IProductWRTBase(), v_NumBndryCoeffs(), and Vmath::Vsub().

◆ v_GenMatrix()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 1221 of file StdTriExp.cpp.

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

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

Referenced by v_CreateStdMatrix().

◆ v_GetBoundaryMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1062 of file StdTriExp.cpp.

1063{
1066 "Expansion not of expected type");
1067 int i;
1068 int cnt;
1069 int nummodes0, nummodes1;
1070 int value;
1071
1072 if (outarray.size() != NumBndryCoeffs())
1073 {
1074 outarray = Array<OneD, unsigned int>(NumBndryCoeffs());
1075 }
1076
1077 nummodes0 = m_base[0]->GetNumModes();
1078 nummodes1 = m_base[1]->GetNumModes();
1079
1080 value = 2 * nummodes1 - 1;
1081 for (i = 0; i < value; i++)
1082 {
1083 outarray[i] = i;
1084 }
1085 cnt = value;
1086
1087 for (i = 0; i < nummodes0 - 2; i++)
1088 {
1089 outarray[cnt++] = value;
1090 value += nummodes1 - 2 - i;
1091 }
1092}

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

◆ v_GetCoords()

void Nektar::StdRegions::StdTriExp::v_GetCoords ( Array< OneD, NekDouble > &  coords_x,
Array< OneD, NekDouble > &  coords_y,
Array< OneD, NekDouble > &  coords_z 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 804 of file StdTriExp.cpp.

807{
808 Array<OneD, const NekDouble> z0 = m_base[0]->GetZ();
809 Array<OneD, const NekDouble> z1 = m_base[1]->GetZ();
810 int nq0 = GetNumPoints(0);
811 int nq1 = GetNumPoints(1);
812 int i, j;
813
814 for (i = 0; i < nq1; ++i)
815 {
816 for (j = 0; j < nq0; ++j)
817 {
818 coords_0[i * nq0 + j] = (1 + z0[j]) * (1 - z1[i]) / 2.0 - 1.0;
819 }
820 Vmath::Fill(nq0, z1[i], &coords_1[0] + i * nq0, 1);
821 }
822}

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

◆ v_GetInteriorMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1032 of file StdTriExp.cpp.

1033{
1036 "Expansion not of a proper type");
1037
1038 int i, j;
1039 int cnt = 0;
1040 int nummodes0, nummodes1;
1041 int startvalue;
1042 if (outarray.size() != GetNcoeffs() - NumBndryCoeffs())
1043 {
1044 outarray = Array<OneD, unsigned int>(GetNcoeffs() - NumBndryCoeffs());
1045 }
1046
1047 nummodes0 = m_base[0]->GetNumModes();
1048 nummodes1 = m_base[1]->GetNumModes();
1049
1050 startvalue = 2 * nummodes1;
1051
1052 for (i = 0; i < nummodes0 - 2; i++)
1053 {
1054 for (j = 0; j < nummodes1 - 3 - i; j++)
1055 {
1056 outarray[cnt++] = startvalue + j;
1057 }
1058 startvalue += nummodes1 - 2 - i;
1059 }
1060}
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_FwdTransBndConstrained().

◆ v_GetNtraces()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 722 of file StdTriExp.cpp.

723{
724 return 3;
725}

◆ v_GetNverts()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 717 of file StdTriExp.cpp.

718{
719 return 3;
720}

◆ v_GetSimplexEquiSpacedConnectivity()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1541 of file StdTriExp.cpp.

1543{
1544 int np1 = m_base[0]->GetNumPoints();
1545 int np2 = m_base[1]->GetNumPoints();
1546 int np = max(np1, np2);
1547
1548 conn = Array<OneD, int>(3 * (np - 1) * (np - 1));
1549
1550 int row = 0;
1551 int rowp1 = 0;
1552 int cnt = 0;
1553 for (int i = 0; i < np - 1; ++i)
1554 {
1555 rowp1 += np - i;
1556 for (int j = 0; j < np - i - 2; ++j)
1557 {
1558 conn[cnt++] = row + j;
1559 conn[cnt++] = row + j + 1;
1560 conn[cnt++] = rowp1 + j;
1561
1562 conn[cnt++] = rowp1 + j + 1;
1563 conn[cnt++] = rowp1 + j;
1564 conn[cnt++] = row + j + 1;
1565 }
1566
1567 conn[cnt++] = row + np - i - 2;
1568 conn[cnt++] = row + np - i - 1;
1569 conn[cnt++] = rowp1 + np - i - 2;
1570
1571 row += np - i;
1572 }
1573}

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

◆ v_GetTraceBasisKey()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 830 of file StdTriExp.cpp.

832{
833 ASSERTL2(i >= 0 && i <= 2, "edge id is out of range");
834
835 // Get basiskey (0 or 1) according to edge id i
836 int dir = (i != 0);
837
838 switch (m_base[dir]->GetBasisType())
839 {
842 {
843 switch (m_base[dir]->GetPointsType())
844 {
846 {
847 return m_base[dir]->GetBasisKey();
848 }
849 break;
850 default:
851 {
853 "Unexpected points distribution " +
855 [m_base[dir]->GetPointsType()] +
856 " in StdTriExp::v_GetTraceBasisKey");
857 }
858 }
859 }
860 break;
862 {
863 switch (m_base[dir]->GetPointsType())
864 {
865 case LibUtilities::eGaussRadauMAlpha1Beta0:
866 {
867 LibUtilities::PointsKey pkey(
868 m_base[dir]
869 ->GetBasisKey()
870 .GetPointsKey()
871 .GetNumPoints() +
872 1,
874 return LibUtilities::BasisKey(LibUtilities::eModified_A,
875 m_base[dir]->GetNumModes(),
876 pkey);
877 }
878 break;
879 // Currently this does not increase the points by
880 // 1 since when using this quadrature we are
881 // presuming it is already been increased by one
882 // when comopared to
883 // GaussRadauMAlpha1Beta0. Currently used in the
884 // GJP option
885 //
886 // Note have put down it back to numpoints +1 to
887 // test for use on tri faces and GJP.
889 {
890 LibUtilities::PointsKey pkey(
891 m_base[dir]
892 ->GetBasisKey()
893 .GetPointsKey()
894 .GetNumPoints() +
895 1,
897 return LibUtilities::BasisKey(LibUtilities::eModified_A,
898 m_base[dir]->GetNumModes(),
899 pkey);
900 }
901 break;
903 {
904 LibUtilities::PointsKey pkey(
905 m_base[dir]
906 ->GetBasisKey()
907 .GetPointsKey()
908 .GetNumPoints() +
909 1,
911 return LibUtilities::BasisKey(LibUtilities::eModified_A,
912 m_base[dir]->GetNumModes(),
913 pkey);
914 }
915 break;
916 default:
917 {
919 "Unexpected points distribution " +
921 [m_base[dir]->GetPointsType()] +
922 " in StdTriExp::v_GetTraceBasisKey");
923 }
924 }
925 }
926 break;
928 {
929 switch (m_base[dir]->GetPointsType())
930 {
931 case LibUtilities::eGaussRadauMAlpha1Beta0:
932 {
933 LibUtilities::PointsKey pkey(
934 m_base[dir]
935 ->GetBasisKey()
936 .GetPointsKey()
937 .GetNumPoints() +
938 1,
940 return LibUtilities::BasisKey(LibUtilities::eGLL_Lagrange,
941 m_base[dir]->GetNumModes(),
942 pkey);
943 }
944 break;
945 default:
946 {
948 "Unexpected points distribution " +
950 [m_base[dir]->GetPointsType()] +
951 " in StdTriExp::v_GetTraceBasisKey");
952 }
953 }
954 }
955 break;
956 default:
957 {
959 "Information not available to set edge key");
960 }
961 }
963}
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:205
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.
const std::string kPointsTypeStr[]
Definition: Foundations.hpp:52
@ eGaussRadauMLegendre
1D Gauss-Radau-Legendre quadrature points, pinned at x=-1
Definition: PointsType.h:47
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
@ ePolyEvenlySpaced
1D Evenly-spaced points using Lagrange polynomial
Definition: PointsType.h:73
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:56

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

◆ v_GetTraceCoeffMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion2D.

Definition at line 1094 of file StdTriExp.cpp.

1096{
1097
1100 "Mapping not defined for this type of basis");
1101
1102 ASSERTL1(eid < 3, "eid must be between 0 and 2");
1103
1104 int i;
1105 int order0 = m_base[0]->GetNumModes();
1106 int order1 = m_base[1]->GetNumModes();
1107 int numModes = (eid == 0) ? order0 : order1;
1108
1109 if (maparray.size() != numModes)
1110 {
1111 maparray = Array<OneD, unsigned int>(numModes);
1112 }
1113
1114 switch (eid)
1115 {
1116 case 0:
1117 {
1118 int cnt = 0;
1119 for (i = 0; i < numModes; cnt += order1 - i, ++i)
1120 {
1121 maparray[i] = cnt;
1122 }
1123 break;
1124 }
1125 case 1:
1126 {
1127 maparray[0] = order1;
1128 maparray[1] = 1;
1129 for (i = 2; i < numModes; i++)
1130 {
1131 maparray[i] = order1 - 1 + i;
1132 }
1133 break;
1134 }
1135 case 2:
1136 {
1137 for (i = 0; i < numModes; i++)
1138 {
1139 maparray[i] = i;
1140 }
1141 break;
1142 }
1143 default:
1144 ASSERTL0(false, "eid must be between 0 and 2");
1145 break;
1146 }
1147}

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

◆ v_GetTraceInteriorToElementMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1149 of file StdTriExp.cpp.

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

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

◆ v_GetTraceIntNcoeffs()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 766 of file StdTriExp.cpp.

767{
768 ASSERTL2(i >= 0 && i <= 2, "edge id is out of range");
769
770 if (i == 0)
771 {
772 return GetBasisNumModes(0) - 2;
773 }
774 else
775 {
776 return GetBasisNumModes(1) - 2;
777 }
778}
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:169

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

◆ v_GetTraceNcoeffs()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 752 of file StdTriExp.cpp.

753{
754 ASSERTL2(i >= 0 && i <= 2, "edge id is out of range");
755
756 if (i == 0)
757 {
758 return GetBasisNumModes(0);
759 }
760 else
761 {
762 return GetBasisNumModes(1);
763 }
764}

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

◆ v_GetTraceNumPoints()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 780 of file StdTriExp.cpp.

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

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

◆ v_GetVertexMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 969 of file StdTriExp.cpp.

970{
973 "Mapping not defined for this type of basis");
974
975 int localDOF = 0;
976 if (useCoeffPacking == true)
977 {
978 switch (localVertexId)
979 {
980 case 0:
981 {
982 localDOF = 0;
983 break;
984 }
985 case 1:
986 {
987 localDOF = 1;
988 break;
989 }
990 case 2:
991 {
992 localDOF = m_base[1]->GetNumModes();
993 break;
994 }
995 default:
996 {
997 ASSERTL0(false, "eid must be between 0 and 2");
998 break;
999 }
1000 }
1001 }
1002 else // follow book format for vertex indexing.
1003 {
1004 switch (localVertexId)
1005 {
1006 case 0:
1007 {
1008 localDOF = 0;
1009 break;
1010 }
1011 case 1:
1012 {
1013 localDOF = m_base[1]->GetNumModes();
1014 break;
1015 }
1016 case 2:
1017 {
1018 localDOF = 1;
1019 break;
1020 }
1021 default:
1022 {
1023 ASSERTL0(false, "eid must be between 0 and 2");
1024 break;
1025 }
1026 }
1027 }
1028
1029 return localDOF;
1030}

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

◆ v_HelmholtzMatrixOp()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 1334 of file StdTriExp.cpp.

1337{
1338 StdTriExp::v_HelmholtzMatrixOp_MatFree(inarray, outarray, mkey);
1339}
void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override

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

◆ v_Integral()

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

Integrates the specified function over the domain.

See also
StdRegions::StdExpansion::Integral.

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 58 of file StdTriExp.cpp.

59{
60 int i;
61 int nquad1 = m_base[1]->GetNumPoints();
62 Array<OneD, NekDouble> w1_tmp(nquad1);
63
64 Array<OneD, const NekDouble> w0 = m_base[0]->GetW();
65 Array<OneD, const NekDouble> w1 = m_base[1]->GetW();
66 Array<OneD, const NekDouble> z1 = m_base[1]->GetZ();
67
68 switch (m_base[1]->GetPointsType())
69 {
70 case LibUtilities::eGaussRadauMAlpha1Beta0: // (0,1) Jacobi Inner
71 // product
72 {
73 Vmath::Smul(nquad1, 0.5, w1, 1, w1_tmp, 1);
74 break;
75 }
76 default:
77 {
78 // include jacobian factor on whatever coordinates are defined.
79 for (i = 0; i < nquad1; ++i)
80 {
81 w1_tmp[i] = 0.5 * (1 - z1[i]) * w1[i];
82 }
83 break;
84 }
85 }
86
87 return StdExpansion2D::Integral(inarray, w0, w1_tmp);
88}
NekDouble Integral(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &w0, const Array< OneD, const NekDouble > &w1)

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

◆ v_IProductWRTBase()

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

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

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

where

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

which can be implemented as

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

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

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

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

Implements Nektar::StdRegions::StdExpansion.

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

Definition at line 420 of file StdTriExp.cpp.

422{
423 StdTriExp::v_IProductWRTBase_SumFac(inarray, outarray);
424}
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
Definition: StdTriExp.cpp:426

References v_IProductWRTBase_SumFac().

Referenced by v_FwdTrans(), and v_FwdTransBndConstrained().

◆ v_IProductWRTBase_SumFac()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 426 of file StdTriExp.cpp.

429{
430 int nquad0 = m_base[0]->GetNumPoints();
431 int nquad1 = m_base[1]->GetNumPoints();
432 int order0 = m_base[0]->GetNumModes();
433
434 if (multiplybyweights)
435 {
436 Array<OneD, NekDouble> tmp(nquad0 * nquad1 + nquad1 * order0);
437 Array<OneD, NekDouble> wsp(tmp + nquad0 * nquad1);
438
439 // multiply by integration constants
440 MultiplyByQuadratureMetric(inarray, tmp);
442 m_base[1]->GetBdata(), tmp, outarray, wsp);
443 }
444 else
445 {
446 Array<OneD, NekDouble> wsp(nquad1 * order0);
448 m_base[1]->GetBdata(), inarray, outarray,
449 wsp);
450 }
451}
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:723

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

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

◆ v_IProductWRTBase_SumFacKernel()

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

Implements Nektar::StdRegions::StdExpansion2D.

Definition at line 453 of file StdTriExp.cpp.

460{
461 int i;
462 int mode;
463 int nquad0 = m_base[0]->GetNumPoints();
464 int nquad1 = m_base[1]->GetNumPoints();
465 int nmodes0 = m_base[0]->GetNumModes();
466 int nmodes1 = m_base[1]->GetNumModes();
467
468 ASSERTL1(wsp.size() >= nquad1 * nmodes0,
469 "Workspace size is not sufficient");
470
471 Blas::Dgemm('T', 'N', nquad1, nmodes0, nquad0, 1.0, inarray.get(), nquad0,
472 base0.get(), nquad0, 0.0, wsp.get(), nquad1);
473
474 // Inner product with respect to 'b' direction
475 for (mode = i = 0; i < nmodes0; ++i)
476 {
477 Blas::Dgemv('T', nquad1, nmodes1 - i, 1.0, base1.get() + mode * nquad1,
478 nquad1, wsp.get() + i * nquad1, 1, 0.0,
479 outarray.get() + mode, 1);
480 mode += nmodes1 - i;
481 }
482
483 // fix for modified basis by splitting top vertex mode
485 {
486 outarray[1] +=
487 Blas::Ddot(nquad1, base1.get() + nquad1, 1, wsp.get() + nquad1, 1);
488 }
489}
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:163

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

◆ v_IProductWRTDerivBase()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 491 of file StdTriExp.cpp.

494{
495 StdTriExp::v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
496}
void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdTriExp.cpp:498

References v_IProductWRTDerivBase_SumFac().

◆ v_IProductWRTDerivBase_SumFac()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 498 of file StdTriExp.cpp.

501{
502 int i;
503 int nquad0 = m_base[0]->GetNumPoints();
504 int nquad1 = m_base[1]->GetNumPoints();
505 int nqtot = nquad0 * nquad1;
506 int nmodes0 = m_base[0]->GetNumModes();
507 int wspsize = max(max(nqtot, m_ncoeffs), nquad1 * nmodes0);
508
509 Array<OneD, NekDouble> gfac0(2 * wspsize);
510 Array<OneD, NekDouble> tmp0(gfac0 + wspsize);
511
512 const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
513
514 // set up geometric factor: 2/(1-z1)
515 for (i = 0; i < nquad1; ++i)
516 {
517 gfac0[i] = 2.0 / (1 - z1[i]);
518 }
519
520 for (i = 0; i < nquad1; ++i)
521 {
522 Vmath::Smul(nquad0, gfac0[i], &inarray[0] + i * nquad0, 1,
523 &tmp0[0] + i * nquad0, 1);
524 }
525
526 MultiplyByQuadratureMetric(tmp0, tmp0);
527
528 switch (dir)
529 {
530 case 0:
531 {
532 IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
533 m_base[1]->GetBdata(), tmp0, outarray,
534 gfac0);
535 break;
536 }
537 case 1:
538 {
539 Array<OneD, NekDouble> tmp3(m_ncoeffs);
540 const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
541
542 for (i = 0; i < nquad0; ++i)
543 {
544 gfac0[i] = 0.5 * (1 + z0[i]);
545 }
546
547 for (i = 0; i < nquad1; ++i)
548 {
549 Vmath::Vmul(nquad0, &gfac0[0], 1, &tmp0[0] + i * nquad0, 1,
550 &tmp0[0] + i * nquad0, 1);
551 }
552
553 IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
554 m_base[1]->GetBdata(), tmp0, tmp3,
555 gfac0);
556
557 MultiplyByQuadratureMetric(inarray, tmp0);
559 m_base[1]->GetDbdata(), tmp0, outarray,
560 gfac0);
561 Vmath::Vadd(m_ncoeffs, &tmp3[0], 1, &outarray[0], 1, &outarray[0],
562 1);
563 break;
564 }
565 default:
566 {
567 ASSERTL1(false, "input dir is out of range");
568 break;
569 }
570 }
571}
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.hpp:180

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

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

◆ v_IsBoundaryInteriorExpansion()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 824 of file StdTriExp.cpp.

825{
826 return m_base[0]->GetBasisType() == LibUtilities::eModified_A &&
827 m_base[1]->GetBasisType() == LibUtilities::eModified_B;
828}

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

◆ v_LaplacianMatrixOp() [1/2]

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 1311 of file StdTriExp.cpp.

1314{
1315 StdTriExp::v_LaplacianMatrixOp_MatFree(inarray, outarray, mkey);
1316}
void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override

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

◆ v_LaplacianMatrixOp() [2/2]

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 1318 of file StdTriExp.cpp.

1322{
1323 StdExpansion::LaplacianMatrixOp_MatFree(k1, k2, inarray, outarray, mkey);
1324}
void LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)

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

◆ v_LocCollapsedToLocCoord()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 596 of file StdTriExp.cpp.

598{
599 xi[0] = (1.0 + eta[0]) * (1.0 - eta[1]) * 0.5 - 1.0;
600 xi[1] = eta[1];
601}

◆ v_LocCoordToLocCollapsed()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 577 of file StdTriExp.cpp.

579{
580 NekDouble d1 = 1. - xi[1];
581 if (fabs(d1) < NekConstants::kNekZeroTol)
582 {
583 if (d1 >= 0.)
584 {
586 }
587 else
588 {
590 }
591 }
592 eta[0] = 2. * (1. + xi[0]) / d1 - 1.0;
593 eta[1] = xi[1];
594}
static const NekDouble kNekZeroTol

References Nektar::NekConstants::kNekZeroTol.

◆ v_MassMatrixOp()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 1304 of file StdTriExp.cpp.

1307{
1308 StdExpansion::MassMatrixOp_MatFree(inarray, outarray, mkey);
1309}
void MassMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)

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

◆ v_MultiplyByStdQuadratureMetric()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1501 of file StdTriExp.cpp.

1504{
1505 int i;
1506 int nquad0 = m_base[0]->GetNumPoints();
1507 int nquad1 = m_base[1]->GetNumPoints();
1508
1509 const Array<OneD, const NekDouble> &w0 = m_base[0]->GetW();
1510 const Array<OneD, const NekDouble> &w1 = m_base[1]->GetW();
1511 const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
1512
1513 // multiply by integration constants
1514 for (i = 0; i < nquad1; ++i)
1515 {
1516 Vmath::Vmul(nquad0, inarray.get() + i * nquad0, 1, w0.get(), 1,
1517 outarray.get() + i * nquad0, 1);
1518 }
1519
1520 switch (m_base[1]->GetPointsType())
1521 {
1522 // (1,0) Jacobi Inner product
1523 case LibUtilities::eGaussRadauMAlpha1Beta0:
1524 for (i = 0; i < nquad1; ++i)
1525 {
1526 Blas::Dscal(nquad0, 0.5 * w1[i], outarray.get() + i * nquad0,
1527 1);
1528 }
1529 break;
1530 // Legendre inner product
1531 default:
1532 for (i = 0; i < nquad1; ++i)
1533 {
1534 Blas::Dscal(nquad0, 0.5 * (1 - z1[i]) * w1[i],
1535 outarray.get() + i * nquad0, 1);
1536 }
1537 break;
1538 }
1539}
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition: Blas.hpp:149

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

◆ v_NumBndryCoeffs()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 732 of file StdTriExp.cpp.

733{
735 "BasisType is not a boundary interior form");
737 "BasisType is not a boundary interior form");
738
739 return 3 + (GetBasisNumModes(0) - 2) + 2 * (GetBasisNumModes(1) - 2);
740}

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

Referenced by v_FwdTransBndConstrained().

◆ v_NumDGBndryCoeffs()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 742 of file StdTriExp.cpp.

743{
745 "BasisType is not a boundary interior form");
747 "BasisType is not a boundary interior form");
748
749 return GetBasisNumModes(0) + 2 * GetBasisNumModes(1);
750}

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

◆ v_PhysDeriv() [1/2]

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

Calculate the derivative of the physical points.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 106 of file StdTriExp.cpp.

110{
111 int i;
112 int nquad0 = m_base[0]->GetNumPoints();
113 int nquad1 = m_base[1]->GetNumPoints();
114 Array<OneD, NekDouble> wsp(std::max(nquad0, nquad1));
115
116 const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
117 const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
118
119 // set up geometric factor: 2/(1-z1)
120 Vmath::Sadd(nquad1, -1.0, z1, 1, wsp, 1);
121 Vmath::Sdiv(nquad1, -2.0, wsp, 1, wsp, 1);
122
123 if (out_d0.size() > 0)
124 {
125 PhysTensorDeriv(inarray, out_d0, out_d1);
126
127 for (i = 0; i < nquad1; ++i)
128 {
129 Blas::Dscal(nquad0, wsp[i], &out_d0[0] + i * nquad0, 1);
130 }
131
132 // if no d1 required do not need to calculate both deriv
133 if (out_d1.size() > 0)
134 {
135 // set up geometric factor: (1+z0)/2
136 Vmath::Sadd(nquad0, 1.0, z0, 1, wsp, 1);
137 Vmath::Smul(nquad0, 0.5, wsp, 1, wsp, 1);
138
139 for (i = 0; i < nquad1; ++i)
140 {
141 Vmath::Vvtvp(nquad0, &wsp[0], 1, &out_d0[0] + i * nquad0, 1,
142 &out_d1[0] + i * nquad0, 1,
143 &out_d1[0] + i * nquad0, 1);
144 }
145 }
146 }
147 else if (out_d1.size() > 0)
148 {
149 Array<OneD, NekDouble> diff0(nquad0 * nquad1);
150 PhysTensorDeriv(inarray, diff0, out_d1);
151
152 for (i = 0; i < nquad1; ++i)
153 {
154 Blas::Dscal(nquad0, wsp[i], &diff0[0] + i * nquad0, 1);
155 }
156
157 Vmath::Sadd(nquad0, 1.0, z0, 1, wsp, 1);
158 Vmath::Smul(nquad0, 0.5, wsp, 1, wsp, 1);
159
160 for (i = 0; i < nquad1; ++i)
161 {
162 Vmath::Vvtvp(nquad0, &wsp[0], 1, &diff0[0] + i * nquad0, 1,
163 &out_d1[0] + i * nquad0, 1, &out_d1[0] + i * nquad0,
164 1);
165 }
166 }
167}
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.hpp:366
void Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/x.
Definition: Vmath.hpp:154
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.hpp:194

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

Referenced by v_PhysDeriv(), and v_StdPhysDeriv().

◆ v_PhysDeriv() [2/2]

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

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

See also
StdRegions::StdExpansion::PhysDeriv

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 169 of file StdTriExp.cpp.

172{
173 switch (dir)
174 {
175 case 0:
176 {
177 v_PhysDeriv(inarray, outarray, NullNekDouble1DArray);
178 break;
179 }
180 case 1:
181 {
182 v_PhysDeriv(inarray, NullNekDouble1DArray, outarray);
183 break;
184 }
185 default:
186 {
187 ASSERTL1(false, "input dir is out of range");
188 break;
189 }
190 }
191}
void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray) override
Calculate the derivative of the physical points.
Definition: StdTriExp.cpp:106
static Array< OneD, NekDouble > NullNekDouble1DArray

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

◆ v_PhysEvaluate()

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

Reimplemented from Nektar::StdRegions::StdExpansion2D.

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

Definition at line 672 of file StdTriExp.cpp.

675{
676 // Collapse coordinates
677 Array<OneD, NekDouble> coll(2, 0.0);
678 LocCoordToLocCollapsed(coord, coll);
679
680 // If near singularity do the old interpolation matrix method
681 if ((1 - coll[1]) < 1e-5)
682 {
683 int totPoints = GetTotPoints();
684 Array<OneD, NekDouble> EphysDeriv0(totPoints), EphysDeriv1(totPoints);
685 PhysDeriv(inarray, EphysDeriv0, EphysDeriv1);
686
687 Array<OneD, DNekMatSharedPtr> I(2);
688 I[0] = GetBase()[0]->GetI(coll);
689 I[1] = GetBase()[1]->GetI(coll + 1);
690
691 firstOrderDerivs[0] = PhysEvaluate(I, EphysDeriv0);
692 firstOrderDerivs[1] = PhysEvaluate(I, EphysDeriv1);
693 return PhysEvaluate(I, inarray);
694 }
695
696 // set up geometric factor: 2.0/(1.0-z1)
697 NekDouble fac0 = 2 / (1 - coll[1]);
698
699 NekDouble val = BaryTensorDeriv(coll, inarray, firstOrderDerivs);
700
701 // Copy d0 into temp for d1
702 NekDouble temp;
703 temp = firstOrderDerivs[0];
704
705 // Multiply by geometric factor
706 firstOrderDerivs[0] = firstOrderDerivs[0] * fac0;
707
708 // set up geometric factor: (1+z0)/(1-z1)
709 NekDouble fac1 = fac0 * (coll[0] + 1) / 2;
710
711 // Multiply out_d0 by geometric factor and add to out_d1
712 firstOrderDerivs[1] += fac1 * temp;
713
714 return val;
715}
NekDouble BaryTensorDeriv(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:134
const Array< OneD, const LibUtilities::BasisSharedPtr > & GetBase() const
This function gets the shared point to basis.
Definition: StdExpansion.h:100
NekDouble PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
This function evaluates the expansion at a single (arbitrary) point of the domain.
Definition: StdExpansion.h:919
void PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Definition: StdExpansion.h:849

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

◆ v_PhysEvaluateBasis()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 648 of file StdTriExp.cpp.

650{
651 Array<OneD, NekDouble> coll(2);
652 LocCoordToLocCollapsed(coords, coll);
653
654 // From mode we need to determine mode0 and mode1 in the (p,q)
655 // direction. mode1 can be directly inferred from mode.
656 const int nm1 = m_base[1]->GetNumModes();
657 const double c = 1 + 2 * nm1;
658 const int mode0 = floor(0.5 * (c - sqrt(c * c - 8 * mode)));
659
660 if (mode == 1 && m_base[0]->GetBasisType() == LibUtilities::eModified_A)
661 {
662 // Account for collapsed vertex.
663 return StdExpansion::BaryEvaluateBasis<1>(coll[1], 1);
664 }
665 else
666 {
667 return StdExpansion::BaryEvaluateBasis<0>(coll[0], mode0) *
668 StdExpansion::BaryEvaluateBasis<1>(coll[1], mode);
669 }
670}
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

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

◆ v_ReduceOrderCoeffs()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp.

Definition at line 1443 of file StdTriExp.cpp.

1446{
1447 int n_coeffs = inarray.size();
1448 int nquad0 = m_base[0]->GetNumPoints();
1449 int nquad1 = m_base[1]->GetNumPoints();
1450 Array<OneD, NekDouble> coeff(n_coeffs);
1451 Array<OneD, NekDouble> coeff_tmp(n_coeffs, 0.0);
1452 Array<OneD, NekDouble> tmp;
1453 Array<OneD, NekDouble> tmp2;
1454 int nqtot = nquad0 * nquad1;
1455 Array<OneD, NekDouble> phys_tmp(nqtot, 0.0);
1456
1457 int nmodes0 = m_base[0]->GetNumModes();
1458 int nmodes1 = m_base[1]->GetNumModes();
1459 int numMin2 = nmodes0;
1460 int i;
1461
1462 const LibUtilities::PointsKey Pkey0(nmodes0,
1464 const LibUtilities::PointsKey Pkey1(nmodes1,
1466
1467 LibUtilities::BasisKey b0(m_base[0]->GetBasisType(), nmodes0, Pkey0);
1468 LibUtilities::BasisKey b1(m_base[1]->GetBasisType(), nmodes1, Pkey1);
1469
1470 LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A, nmodes0, Pkey0);
1471 LibUtilities::BasisKey bortho1(LibUtilities::eOrtho_B, nmodes1, Pkey1);
1472
1473 StdRegions::StdTriExpSharedPtr m_OrthoTriExp;
1475
1478 bortho0, bortho1);
1479
1480 m_TriExp->BwdTrans(inarray, phys_tmp);
1481 m_OrthoTriExp->FwdTrans(phys_tmp, coeff);
1482
1483 for (i = 0; i < n_coeffs; i++)
1484 {
1485 if (i == numMin)
1486 {
1487 coeff[i] = 0.0;
1488 numMin += numMin2 - 1;
1489 numMin2 -= 1.0;
1490 }
1491 }
1492
1493 m_OrthoTriExp->BwdTrans(coeff, phys_tmp);
1494 m_TriExp->FwdTrans(phys_tmp, outarray);
1495}
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:42
std::shared_ptr< StdTriExp > StdTriExpSharedPtr
Definition: StdTriExp.h:219

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

◆ v_StdPhysDeriv() [1/2]

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::NodalTriExp.

Definition at line 193 of file StdTriExp.cpp.

197{
198 StdTriExp::v_PhysDeriv(inarray, out_d0, out_d1);
199}

References v_PhysDeriv().

◆ v_StdPhysDeriv() [2/2]

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 201 of file StdTriExp.cpp.

204{
205 StdTriExp::v_PhysDeriv(dir, inarray, outarray);
206}

References v_PhysDeriv().

◆ v_SVVLaplacianFilter()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp.

Definition at line 1341 of file StdTriExp.cpp.

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

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

◆ v_WeakDerivMatrixOp()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 1326 of file StdTriExp.cpp.

1330{
1331 StdExpansion::WeakDerivMatrixOp_MatFree(i, inarray, outarray, mkey);
1332}
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().