Nektar++
Public Member Functions | Protected Member Functions | Private Member Functions | List of all members
Nektar::StdRegions::StdExpansion2D Class Referenceabstract

#include <StdExpansion2D.h>

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

Public Member Functions

 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, bool UseGLL=false) 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 NodalToModal (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
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...
 
void PhysInterp (std::shared_ptr< StdExpansion > fromExp, const Array< OneD, const NekDouble > &fromData, Array< OneD, NekDouble > &toData)
 interpolate from one set of quadrature points available from FromExp to the set of quadrature points in the current expansion. If the points are the same this routine will just copy the data 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_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_PhysEvaluateInterp (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals) 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
 
void v_PhysInterp (std::shared_ptr< StdExpansion > fromExp, const Array< OneD, const NekDouble > &fromData, Array< OneD, NekDouble > &toData) 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 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)
 

Private Member Functions

int v_GetShapeDimension () const final
 

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 44 of file StdExpansion2D.h.

Constructor & Destructor Documentation

◆ StdExpansion2D() [1/3]

Nektar::StdRegions::StdExpansion2D::StdExpansion2D ( int  numcoeffs,
const LibUtilities::BasisKey Ba,
const LibUtilities::BasisKey Bb 
)

Definition at line 47 of file StdExpansion2D.cpp.

51{
52}

◆ StdExpansion2D() [2/3]

Nektar::StdRegions::StdExpansion2D::StdExpansion2D ( )
default

◆ StdExpansion2D() [3/3]

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

◆ ~StdExpansion2D()

Nektar::StdRegions::StdExpansion2D::~StdExpansion2D ( )
overridedefault

Member Function Documentation

◆ BaryTensorDeriv()

NekDouble Nektar::StdRegions::StdExpansion2D::BaryTensorDeriv ( const Array< OneD, NekDouble > &  coord,
const Array< OneD, const NekDouble > &  inarray,
std::array< NekDouble, 3 > &  firstOrderDerivs 
)
inline

Definition at line 99 of file StdExpansion2D.h.

103 {
104 const int nq0 = m_base[0]->GetNumPoints();
105 const int nq1 = m_base[1]->GetNumPoints();
106
107 const NekDouble *ptr = &inarray[0];
108 Array<OneD, NekDouble> deriv0(nq1, 0.0);
109 Array<OneD, NekDouble> phys0(nq1, 0.0);
110
111 for (int j = 0; j < nq1; ++j, ptr += nq0)
112 {
113 phys0[j] =
114 StdExpansion::BaryEvaluate<0, true>(coord[0], ptr, deriv0[j]);
115 }
116 firstOrderDerivs[0] =
117 StdExpansion::BaryEvaluate<1, false>(coord[1], &deriv0[0]);
118
119 return StdExpansion::BaryEvaluate<1, true>(coord[1], &phys0[0],
120 firstOrderDerivs[1]);
121 }
Array< OneD, LibUtilities::BasisSharedPtr > m_base
double NekDouble

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

Referenced by Nektar::StdRegions::StdQuadExp::v_PhysEvalFirstDeriv(), and Nektar::StdRegions::StdTriExp::v_PhysEvalFirstDeriv().

◆ BwdTrans_SumFacKernel()

void Nektar::StdRegions::StdExpansion2D::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 
)

Definition at line 181 of file StdExpansion2D.cpp.

187{
188 v_BwdTrans_SumFacKernel(base0, base1, inarray, outarray, wsp,
189 doCheckCollDir0, doCheckCollDir1);
190}
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

References v_BwdTrans_SumFacKernel().

Referenced by Nektar::StdRegions::StdQuadExp::v_BwdTrans_SumFac(), Nektar::StdRegions::StdTriExp::v_BwdTrans_SumFac(), v_HelmholtzMatrixOp_MatFree(), and v_LaplacianMatrixOp_MatFree().

◆ Integral()

NekDouble Nektar::StdRegions::StdExpansion2D::Integral ( const Array< OneD, const NekDouble > &  inarray,
const Array< OneD, const NekDouble > &  w0,
const Array< OneD, const NekDouble > &  w1 
)

Definition at line 154 of file StdExpansion2D.cpp.

157{
158 int i;
159 NekDouble Int = 0.0;
160 int nquad0 = m_base[0]->GetNumPoints();
161 int nquad1 = m_base[1]->GetNumPoints();
162 Array<OneD, NekDouble> tmp(nquad0 * nquad1);
163
164 // multiply by integration constants
165 for (i = 0; i < nquad1; ++i)
166 {
167 Vmath::Vmul(nquad0, &inarray[0] + i * nquad0, 1, w0.data(), 1,
168 &tmp[0] + i * nquad0, 1);
169 }
170
171 for (i = 0; i < nquad0; ++i)
172 {
173 Vmath::Vmul(nquad1, &tmp[0] + i, nquad0, w1.data(), 1, &tmp[0] + i,
174 nquad0);
175 }
176 Int = Vmath::Vsum(nquad0 * nquad1, tmp, 1);
177
178 return Int;
179}
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
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.hpp:608

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

Referenced by Nektar::StdRegions::StdQuadExp::v_Integral(), and Nektar::StdRegions::StdTriExp::v_Integral().

◆ IProductWRTBase_SumFacKernel()

void Nektar::StdRegions::StdExpansion2D::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 
)

◆ PhysTensorDeriv()

void Nektar::StdRegions::StdExpansion2D::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.

This function is independent of the expansion basis and can therefore be defined for all tensor product distribution of quadrature points in a generic manner. The key operations are:

  • \( \frac{d}{d\eta_1} \rightarrow {\bf D^T_0 u } \)
  • \( \frac{d}{d\eta_2} \rightarrow {\bf D_1 u } \)
Parameters
inarrayarray of physical points to be differentiated
outarray_d0the resulting array of derivative in the \(\eta_1\) direction will be stored in outarray_d0 as output of the function
outarray_d1the resulting array of derivative in the \(\eta_2\) direction will be stored in outarray_d1 as output of the function

Recall that: \( \hspace{1cm} \begin{array}{llll} \mbox{Shape} & \mbox{Cartesian coordinate range} & \mbox{Collapsed coord.} & \mbox{Collapsed coordinate definition}\\ \mbox{Quadrilateral} & -1 \leq \xi_1,\xi_2 \leq 1 & -1 \leq \eta_1,\eta_2 \leq 1 & \eta_1 = \xi_1, \eta_2 = \xi_2\\ \mbox{Triangle} & -1 \leq \xi_1,\xi_2; \xi_1+\xi_2 \leq 0 & -1 \leq \eta_1,\eta_2 \leq 1 & \eta_1 = \frac{2(1+\xi_1)}{(1-\xi_2)}-1, \eta_2 = \xi_2 \\ \end{array} \)

Definition at line 57 of file StdExpansion2D.cpp.

60{
61 int nquad0 = m_base[0]->GetNumPoints();
62 int nquad1 = m_base[1]->GetNumPoints();
63
64 if (outarray_d0.size() > 0) // calculate du/dx_0
65 {
66 DNekMatSharedPtr D0 = m_base[0]->GetD();
67 if (inarray.data() == outarray_d0.data())
68 {
69 Array<OneD, NekDouble> wsp(nquad0 * nquad1);
70 Vmath::Vcopy(nquad0 * nquad1, inarray.data(), 1, wsp.data(), 1);
71 Blas::Dgemm('N', 'N', nquad0, nquad1, nquad0, 1.0,
72 &(D0->GetPtr())[0], nquad0, &wsp[0], nquad0, 0.0,
73 &outarray_d0[0], nquad0);
74 }
75 else
76 {
77 Blas::Dgemm('N', 'N', nquad0, nquad1, nquad0, 1.0,
78 &(D0->GetPtr())[0], nquad0, &inarray[0], nquad0, 0.0,
79 &outarray_d0[0], nquad0);
80 }
81 }
82
83 if (outarray_d1.size() > 0) // calculate du/dx_1
84 {
85 DNekMatSharedPtr D1 = m_base[1]->GetD();
86 if (inarray.data() == outarray_d1.data())
87 {
88 Array<OneD, NekDouble> wsp(nquad0 * nquad1);
89 Vmath::Vcopy(nquad0 * nquad1, inarray.data(), 1, wsp.data(), 1);
90 Blas::Dgemm('N', 'T', nquad0, nquad1, nquad1, 1.0, &wsp[0], nquad0,
91 &(D1->GetPtr())[0], nquad1, 0.0, &outarray_d1[0],
92 nquad0);
93 }
94 else
95 {
96 Blas::Dgemm('N', 'T', nquad0, nquad1, nquad1, 1.0, &inarray[0],
97 nquad0, &(D1->GetPtr())[0], nquad1, 0.0,
98 &outarray_d1[0], nquad0);
99 }
100 }
101}
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
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References Blas::Dgemm(), Nektar::StdRegions::StdExpansion::m_base, and Vmath::Vcopy().

Referenced by Nektar::StdRegions::StdQuadExp::v_PhysDeriv(), and Nektar::StdRegions::StdTriExp::v_PhysDeriv().

◆ v_BwdTrans_SumFacKernel()

virtual void Nektar::StdRegions::StdExpansion2D::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 
)
protectedpure virtual

◆ v_GenStdMatBwdDeriv()

void Nektar::StdRegions::StdExpansion2D::v_GenStdMatBwdDeriv ( const int  dir,
DNekMatSharedPtr mat 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 203 of file StdExpansion2D.cpp.

204{
205 ASSERTL1((dir == 0) || (dir == 1), "Invalid direction.");
206
207 int nquad0 = m_base[0]->GetNumPoints();
208 int nquad1 = m_base[1]->GetNumPoints();
209 int nqtot = nquad0 * nquad1;
210 int nmodes0 = m_base[0]->GetNumModes();
211
212 Array<OneD, NekDouble> tmp1(2 * nqtot + m_ncoeffs + nmodes0 * nquad1, 0.0);
213 Array<OneD, NekDouble> tmp3(tmp1 + 2 * nqtot);
214 Array<OneD, NekDouble> tmp4(tmp1 + 2 * nqtot + m_ncoeffs);
215
216 switch (dir)
217 {
218 case 0:
219 for (int i = 0; i < nqtot; i++)
220 {
221 tmp1[i] = 1.0;
222 IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
223 m_base[1]->GetBdata(), tmp1, tmp3,
224 tmp4, false, true);
225 tmp1[i] = 0.0;
226
227 for (int j = 0; j < m_ncoeffs; j++)
228 {
229 (*mat)(j, i) = tmp3[j];
230 }
231 }
232 break;
233 case 1:
234 for (int i = 0; i < nqtot; i++)
235 {
236 tmp1[i] = 1.0;
238 m_base[1]->GetDbdata(), tmp1, tmp3,
239 tmp4, true, false);
240 tmp1[i] = 0.0;
241
242 for (int j = 0; j < m_ncoeffs; j++)
243 {
244 (*mat)(j, i) = tmp3[j];
245 }
246 }
247 break;
248 default:
249 NEKERROR(ErrorUtil::efatal, "Not a 2D expansion.");
250 break;
251 }
252}
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
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)

References ASSERTL1, Nektar::ErrorUtil::efatal, IProductWRTBase_SumFacKernel(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, and NEKERROR.

◆ v_GetElmtTraceToTraceMap()

void Nektar::StdRegions::StdExpansion2D::v_GetElmtTraceToTraceMap ( const unsigned int  eid,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray,
Orientation  edgeOrient,
int  P,
int  Q 
)
overrideprotectedvirtual

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.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 376 of file StdExpansion2D.cpp.

380{
381 unsigned int i;
382
383 int dir;
384 // determine basis direction for edge.
386 {
387 dir = (eid == 0) ? 0 : 1;
388 }
389 else
390 {
391 dir = eid % 2;
392 }
393
394 int numModes = m_base[dir]->GetNumModes();
395
396 // P is the desired length of the map
397 P = (P == -1) ? numModes : P;
398
399 // decalare maparray
400 if (maparray.size() != P)
401 {
402 maparray = Array<OneD, unsigned int>(P);
403 }
404
405 // fill default mapping as increasing index
406 for (i = 0; i < P; ++i)
407 {
408 maparray[i] = i;
409 }
410
411 if (signarray.size() != P)
412 {
413 signarray = Array<OneD, int>(P, 1);
414 }
415 else
416 {
417 std::fill(signarray.data(), signarray.data() + P, 1);
418 }
419
420 // Zero signmap and set maparray to zero if
421 // elemental modes are not as large as trace modes
422 for (i = numModes; i < P; ++i)
423 {
424 signarray[i] = 0.0;
425 maparray[i] = maparray[0];
426 }
427
428 if (edgeOrient == eBackwards)
429 {
430 const LibUtilities::BasisType bType = GetBasisType(dir);
431
432 if ((bType == LibUtilities::eModified_A) ||
433 (bType == LibUtilities::eModified_B))
434 {
435 std::swap(maparray[0], maparray[1]);
436
437 for (i = 3; i < std::min(P, numModes); i += 2)
438 {
439 signarray[i] *= -1;
440 }
441 }
442 else if (bType == LibUtilities::eGLL_Lagrange ||
444 {
445 ASSERTL1(P == numModes, "Different trace space edge dimension "
446 "and element edge dimension not currently "
447 "possible for GLL-Lagrange bases");
448
449 std::reverse(maparray.data(), maparray.data() + P);
450 }
451 else
452 {
453 ASSERTL0(false, "Mapping not defined for this type of basis");
454 }
455 }
456}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:156
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:370
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:49
@ P
Monomial polynomials .
Definition: BasisType.h:62
@ eGauss_Lagrange
Lagrange Polynomials using the Gauss points.
Definition: BasisType.h:57
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:56
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48

References ASSERTL0, ASSERTL1, Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::StdRegions::eBackwards, Nektar::LibUtilities::eGauss_Lagrange, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::LibUtilities::eTriangle, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::m_base, and Nektar::LibUtilities::P.

Referenced by v_GetTraceToElementMap().

◆ v_GetShapeDimension()

int Nektar::StdRegions::StdExpansion2D::v_GetShapeDimension ( ) const
inlinefinalprivatevirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 214 of file StdExpansion2D.h.

215 {
216 return 2;
217 }

◆ v_GetTraceCoeffMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdQuadExp, and Nektar::StdRegions::StdTriExp.

Definition at line 363 of file StdExpansion2D.cpp.

366{
367 ASSERTL0(false,
368 "This method must be defined at the individual shape level");
369}

References ASSERTL0.

Referenced by v_GetTraceToElementMap().

◆ v_GetTraceToElementMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp.

Definition at line 458 of file StdExpansion2D.cpp.

463{
464 Array<OneD, unsigned int> map1, map2;
465 v_GetTraceCoeffMap(eid, map1);
466 v_GetElmtTraceToTraceMap(eid, map2, signarray, edgeOrient, P, Q);
467
468 if (maparray.size() != map2.size())
469 {
470 maparray = Array<OneD, unsigned int>(map2.size());
471 }
472
473 for (int i = 0; i < map2.size(); ++i)
474 {
475 maparray[i] = map1[map2[i]];
476 }
477}
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...

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

◆ v_HelmholtzMatrixOp_MatFree()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 303 of file StdExpansion2D.cpp.

306{
307 if (mkey.GetNVarCoeff() == 0 &&
308 !mkey.ConstFactorExists(StdRegions::eFactorCoeffD00) &&
309 !mkey.ConstFactorExists(StdRegions::eFactorSVVCutoffRatio))
310 {
311 using std::max;
312
313 int nquad0 = m_base[0]->GetNumPoints();
314 int nquad1 = m_base[1]->GetNumPoints();
315 int nqtot = nquad0 * nquad1;
316 int nmodes0 = m_base[0]->GetNumModes();
317 int nmodes1 = m_base[1]->GetNumModes();
318 int wspsize =
319 max(max(max(nqtot, m_ncoeffs), nquad1 * nmodes0), nquad0 * nmodes1);
320 NekDouble lambda = mkey.GetConstFactor(StdRegions::eFactorLambda);
321
322 const Array<OneD, const NekDouble> &base0 = m_base[0]->GetBdata();
323 const Array<OneD, const NekDouble> &base1 = m_base[1]->GetBdata();
324
325 // Allocate temporary storage
326 Array<OneD, NekDouble> wsp0(5 * wspsize); // size wspsize
327 Array<OneD, NekDouble> wsp1(wsp0 + wspsize); // size wspsize
328 Array<OneD, NekDouble> wsp2(wsp0 + 2 * wspsize); // size 3*wspsize
329
330 if (!(m_base[0]->Collocation() && m_base[1]->Collocation()))
331 {
332 // MASS MATRIX OPERATION
333 // The following is being calculated:
334 // wsp0 = B * u_hat = u
335 // wsp1 = W * wsp0
336 // outarray = B^T * wsp1 = B^T * W * B * u_hat = M * u_hat
337 BwdTrans_SumFacKernel(base0, base1, inarray, wsp0, wsp2, true,
338 true);
339 MultiplyByQuadratureMetric(wsp0, wsp1);
340 IProductWRTBase_SumFacKernel(base0, base1, wsp1, outarray, wsp2,
341 true, true);
342
343 LaplacianMatrixOp_MatFree_Kernel(wsp0, wsp1, wsp2);
344 }
345 else
346 {
347 MultiplyByQuadratureMetric(inarray, outarray);
348 LaplacianMatrixOp_MatFree_Kernel(inarray, wsp1, wsp2);
349 }
350
351 // outarray = lambda * outarray + wsp1
352 // = (lambda * M + L ) * u_hat
353 Vmath::Svtvp(m_ncoeffs, lambda, &outarray[0], 1, &wsp1[0], 1,
354 &outarray[0], 1);
355 }
356 else
357 {
359 mkey);
360 }
361}
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 MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:732
void HelmholtzMatrixOp_MatFree_GenericImpl(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 Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Svtvp (scalar times vector plus vector): z = alpha*x + y.
Definition: Vmath.hpp:396

References BwdTrans_SumFacKernel(), Nektar::StdRegions::StdMatrixKey::ConstFactorExists(), Nektar::StdRegions::eFactorCoeffD00, Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdMatrixKey::GetNVarCoeff(), Nektar::StdRegions::StdExpansion::HelmholtzMatrixOp_MatFree_GenericImpl(), IProductWRTBase_SumFacKernel(), Nektar::StdRegions::StdExpansion::LaplacianMatrixOp_MatFree_Kernel(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::StdRegions::StdExpansion::MultiplyByQuadratureMetric(), and Vmath::Svtvp().

Referenced by Nektar::StdRegions::StdQuadExp::v_HelmholtzMatrixOp(), and Nektar::StdRegions::StdTriExp::v_HelmholtzMatrixOp().

◆ v_IProductWRTBase_SumFacKernel()

virtual void Nektar::StdRegions::StdExpansion2D::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 
)
protectedpure virtual

◆ v_LaplacianMatrixOp_MatFree()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 254 of file StdExpansion2D.cpp.

257{
258 if (mkey.GetNVarCoeff() == 0 &&
259 !mkey.ConstFactorExists(StdRegions::eFactorCoeffD00) &&
260 !mkey.ConstFactorExists(StdRegions::eFactorSVVCutoffRatio))
261 {
262 using std::max;
263
264 // This implementation is only valid when there are no
265 // coefficients associated to the Laplacian operator
266 int nquad0 = m_base[0]->GetNumPoints();
267 int nquad1 = m_base[1]->GetNumPoints();
268 int nqtot = nquad0 * nquad1;
269 int nmodes0 = m_base[0]->GetNumModes();
270 int nmodes1 = m_base[1]->GetNumModes();
271 int wspsize =
272 max(max(max(nqtot, m_ncoeffs), nquad1 * nmodes0), nquad0 * nmodes1);
273
274 const Array<OneD, const NekDouble> &base0 = m_base[0]->GetBdata();
275 const Array<OneD, const NekDouble> &base1 = m_base[1]->GetBdata();
276
277 // Allocate temporary storage
278 Array<OneD, NekDouble> wsp0(4 * wspsize); // size wspsize
279 Array<OneD, NekDouble> wsp1(wsp0 + wspsize); // size 3*wspsize
280
281 if (!(m_base[0]->Collocation() && m_base[1]->Collocation()))
282 {
283 // LAPLACIAN MATRIX OPERATION
284 // wsp0 = u = B * u_hat
285 // wsp1 = du_dxi1 = D_xi1 * wsp0 = D_xi1 * u
286 // wsp2 = du_dxi2 = D_xi2 * wsp0 = D_xi2 * u
287 BwdTrans_SumFacKernel(base0, base1, inarray, wsp0, wsp1, true,
288 true);
289 LaplacianMatrixOp_MatFree_Kernel(wsp0, outarray, wsp1);
290 }
291 else
292 {
293 LaplacianMatrixOp_MatFree_Kernel(inarray, outarray, wsp1);
294 }
295 }
296 else
297 {
299 mkey);
300 }
301}
void LaplacianMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)

References BwdTrans_SumFacKernel(), Nektar::StdRegions::StdMatrixKey::ConstFactorExists(), Nektar::StdRegions::eFactorCoeffD00, Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::StdMatrixKey::GetNVarCoeff(), Nektar::StdRegions::StdExpansion::LaplacianMatrixOp_MatFree_GenericImpl(), Nektar::StdRegions::StdExpansion::LaplacianMatrixOp_MatFree_Kernel(), Nektar::StdRegions::StdExpansion::m_base, and Nektar::StdRegions::StdExpansion::m_ncoeffs.

Referenced by Nektar::StdRegions::StdQuadExp::v_LaplacianMatrixOp(), and Nektar::StdRegions::StdTriExp::v_LaplacianMatrixOp().

◆ v_PhysEvaluate()

NekDouble Nektar::StdRegions::StdExpansion2D::v_PhysEvaluate ( const Array< OneD, const NekDouble > &  coords,
const Array< OneD, const NekDouble > &  physvals 
)
overrideprotectedvirtual

This function evaluates the expansion at a single (arbitrary) point of the domain.

This function is a wrapper around the virtual function v_PhysEvaluate()

Based on the value of the expansion at the quadrature points, this function calculates the value of the expansion at an arbitrary single points (with coordinates \( \mathbf{x_c}\) given by the pointer coords). This operation, equivalent to

\[ u(\mathbf{x_c}) = \sum_p \phi_p(\mathbf{x_c}) \hat{u}_p \]

is evaluated using Lagrangian interpolants through the quadrature points:

\[ u(\mathbf{x_c}) = \sum_p h_p(\mathbf{x_c}) u_p\]

This function requires that the physical value array \(\mathbf{u}\) (implemented as the attribute #m_phys) is set.

Parameters
coordsthe coordinates of the single point
Returns
returns the value of the expansion at the single point

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 103 of file StdExpansion2D.cpp.

106{
107 ASSERTL2(coords[0] > -1 - NekConstants::kNekZeroTol, "coord[0] < -1");
108 ASSERTL2(coords[0] < 1 + NekConstants::kNekZeroTol, "coord[0] > 1");
109 ASSERTL2(coords[1] > -1 - NekConstants::kNekZeroTol, "coord[1] < -1");
110 ASSERTL2(coords[1] < 1 + NekConstants::kNekZeroTol, "coord[1] > 1");
111
112 Array<OneD, NekDouble> coll(2);
113 LocCoordToLocCollapsed(coords, coll);
114
115 const int nq0 = m_base[0]->GetNumPoints();
116 const int nq1 = m_base[1]->GetNumPoints();
117
118 Array<OneD, NekDouble> wsp(nq1);
119 for (int i = 0; i < nq1; ++i)
120 {
121 wsp[i] = StdExpansion::BaryEvaluate<0>(coll[0], &physvals[0] + i * nq0);
122 }
123
124 return StdExpansion::BaryEvaluate<1>(coll[1], &wsp[0]);
125}
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:265
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
static const NekDouble kNekZeroTol

References ASSERTL2, Nektar::NekConstants::kNekZeroTol, Nektar::StdRegions::StdExpansion::LocCoordToLocCollapsed(), and Nektar::StdRegions::StdExpansion::m_base.

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

◆ v_PhysEvaluateInterp()

NekDouble Nektar::StdRegions::StdExpansion2D::v_PhysEvaluateInterp ( const Array< OneD, DNekMatSharedPtr > &  I,
const Array< OneD, const NekDouble > &  physvals 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 127 of file StdExpansion2D.cpp.

130{
131 NekDouble val;
132 int i;
133 int nq0 = m_base[0]->GetNumPoints();
134 int nq1 = m_base[1]->GetNumPoints();
135 Array<OneD, NekDouble> wsp1(nq1);
136
137 // interpolate first coordinate direction
138 for (i = 0; i < nq1; ++i)
139 {
140 wsp1[i] =
141 Blas::Ddot(nq0, &(I[0]->GetPtr())[0], 1, &physvals[i * nq0], 1);
142 }
143
144 // interpolate in second coordinate direction
145 val = Blas::Ddot(nq1, I[1]->GetPtr(), 1, wsp1, 1);
146
147 return val;
148}
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 Blas::Ddot(), and Nektar::StdRegions::StdExpansion::m_base.

◆ v_PhysInterp()

void Nektar::StdRegions::StdExpansion2D::v_PhysInterp ( std::shared_ptr< StdExpansion fromExp,
const Array< OneD, const NekDouble > &  fromData,
Array< OneD, NekDouble > &  toData 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 479 of file StdExpansion2D.cpp.

482{
483
484 LibUtilities::Interp2D(fromExp->GetBasis(0)->GetPointsKey(),
485 fromExp->GetBasis(1)->GetPointsKey(), fromData,
486 m_base[0]->GetPointsKey(), m_base[1]->GetPointsKey(),
487 toData);
488}
void Interp2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
this function interpolates a 2D function evaluated at the quadrature points of the 2D basis,...
Definition: Interp.cpp:101

References Nektar::LibUtilities::Interp2D(), and Nektar::StdRegions::StdExpansion::m_base.