Nektar++
Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes | List of all members
Nektar::LocalRegions::SegExp Class Reference

#include <SegExp.h>

Inheritance diagram for Nektar::LocalRegions::SegExp:
[legend]

Public Member Functions

 SegExp (const LibUtilities::BasisKey &Ba, const SpatialDomains::Geometry1DSharedPtr &geom)
 Constructor using BasisKey class for quadrature points and order definition. More...
 
 SegExp (const SegExp &S)
 Copy Constructor. More...
 
 ~SegExp () override=default
 
- Public Member Functions inherited from Nektar::StdRegions::StdSegExp
 StdSegExp (const LibUtilities::BasisKey &Ba)
 Constructor using BasisKey class for quadrature points and order definition. More...
 
 StdSegExp (const StdSegExp &T)=default
 
 ~StdSegExp () override=default
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion1D
 StdExpansion1D (int numcoeffs, const LibUtilities::BasisKey &Ba)
 
 StdExpansion1D ()=default
 
 StdExpansion1D (const StdExpansion1D &T)=default
 
 ~StdExpansion1D () override=default
 
void PhysTensorDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Evaluate the derivative \( d/d{\xi_1} \) at the physical quadrature points given by inarray and return in outarray. More...
 
NekDouble BaryTensorDeriv (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
 
NekDouble BaryTensorDeriv (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs, std::array< NekDouble, 6 > &secondOrderDerivs)
 
- 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)
 
- Public Member Functions inherited from Nektar::LocalRegions::Expansion1D
 Expansion1D (SpatialDomains::Geometry1DSharedPtr pGeom)
 
 ~Expansion1D () override=default
 
void AddNormTraceInt (const int dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void AddHDGHelmholtzTraceTerms (const NekDouble tau, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
SpatialDomains::Geometry1DSharedPtr GetGeom1D () const
 
- Public Member Functions inherited from Nektar::LocalRegions::Expansion
 Expansion (SpatialDomains::GeometrySharedPtr pGeom)
 
 Expansion (const Expansion &pSrc)
 
 ~Expansion () override
 
void SetTraceExp (const int traceid, ExpansionSharedPtr &f)
 
ExpansionSharedPtr GetTraceExp (const int traceid)
 
DNekScalMatSharedPtr GetLocMatrix (const LocalRegions::MatrixKey &mkey)
 
void DropLocMatrix (const LocalRegions::MatrixKey &mkey)
 
DNekScalMatSharedPtr GetLocMatrix (const StdRegions::MatrixType mtype, const StdRegions::ConstFactorMap &factors=StdRegions::NullConstFactorMap, const StdRegions::VarCoeffMap &varcoeffs=StdRegions::NullVarCoeffMap)
 
SpatialDomains::GeometrySharedPtr GetGeom () const
 
void Reset ()
 
IndexMapValuesSharedPtr CreateIndexMap (const IndexMapKey &ikey)
 
DNekScalBlkMatSharedPtr CreateStaticCondMatrix (const MatrixKey &mkey)
 
const SpatialDomains::GeomFactorsSharedPtrGetMetricInfo () const
 
DNekMatSharedPtr BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
 
DNekMatSharedPtr BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd)
 
void ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType)
 
void AddEdgeNormBoundaryInt (const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
void AddEdgeNormBoundaryInt (const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
void AddFaceNormBoundaryInt (const int face, const std::shared_ptr< Expansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
void DGDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &coeffs, Array< OneD, NekDouble > &outarray)
 
NekDouble VectorFlux (const Array< OneD, Array< OneD, NekDouble > > &vec)
 
void NormalTraceDerivFactors (Array< OneD, Array< OneD, NekDouble > > &factors, Array< OneD, Array< OneD, NekDouble > > &d0factors, Array< OneD, Array< OneD, NekDouble > > &d1factors)
 
IndexMapValuesSharedPtr GetIndexMap (const IndexMapKey &ikey)
 
void AlignVectorToCollapsedDir (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
ExpansionSharedPtr GetLeftAdjacentElementExp () const
 
ExpansionSharedPtr GetRightAdjacentElementExp () const
 
int GetLeftAdjacentElementTrace () const
 
int GetRightAdjacentElementTrace () const
 
void SetAdjacentElementExp (int traceid, ExpansionSharedPtr &e)
 
StdRegions::Orientation GetTraceOrient (int trace)
 
void SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void DivideByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Divided by the metric jacobi and quadrature weights. More...
 
void GetTraceQFactors (const int trace, Array< OneD, NekDouble > &outarray)
 Extract the metric factors to compute the contravariant fluxes along edge edge and stores them into outarray following the local edge orientation (i.e. anticlockwise convention). More...
 
void GetTracePhysVals (const int trace, const StdRegions::StdExpansionSharedPtr &TraceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=StdRegions::eNoOrientation)
 
void GetTracePhysMap (const int edge, Array< OneD, int > &outarray)
 
void ReOrientTracePhysMap (const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1)
 
const NormalVectorGetTraceNormal (const int id)
 
void ComputeTraceNormal (const int id)
 
const Array< OneD, const NekDouble > & GetPhysNormals (void)
 
void SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
void SetUpPhysNormals (const int trace)
 
void AddRobinMassMatrix (const int traceid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
void TraceNormLen (const int traceid, NekDouble &h, NekDouble &p)
 
void AddRobinTraceContribution (const int traceid, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs)
 
const Array< OneD, const NekDouble > & GetElmtBndNormDirElmtLen (const int nbnd) const
 
void StdDerivBaseOnTraceMat (Array< OneD, DNekMatSharedPtr > &DerivMat)
 

Protected Member Functions

NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray) override
 Integrate the physical point list inarray over region and return the value. More...
 
void v_PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray) override
 Evaluate the derivative \( d/d{\xi_1} \) at the physical quadrature points given by inarray and return in outarray. 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_PhysDeriv_s (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds) override
 Evaluate the derivative along a line: \( d/ds=\frac{spacedim}{||tangent||}d/d{\xi} \). The derivative is calculated performing the product \( du/d{s}=\nabla u \cdot tangent \). More...
 
void v_PhysDeriv_n (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn) override
 Evaluate the derivative normal to a line: \( d/dn=\frac{spacedim}{||normal||}d/d{\xi} \). The derivative is calculated performing the product \( du/d{s}=\nabla u \cdot normal \). More...
 
void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Forward transform from physical quadrature space stored in inarray and evaluate the expansion coefficients and store in outarray. 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
 Inner product of inarray over region with respect to the expansion basis (this)->_Base[0] and return in outarray. More...
 
void v_IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check) override
 Inner product of inarray over region with respect to expansion basis base and return in outarray. More...
 
void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray) override
 
void v_NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray) override
 
NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals) override
 
NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &coord, 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
 
NekDouble v_PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs, std::array< NekDouble, 6 > &secondOrderDerivs) override
 
void v_GetCoord (const Array< OneD, const NekDouble > &Lcoords, Array< OneD, NekDouble > &coords) override
 
void v_GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
 
void v_GetVertexPhysVals (const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray) override
 
void v_GetTracePhysVals (const int edge, const StdRegions::StdExpansionSharedPtr &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient) override
 
void v_GetTracePhysMap (const int vertex, Array< OneD, int > &map) override
 
StdRegions::StdExpansionSharedPtr v_GetStdExp (void) const override
 
StdRegions::StdExpansionSharedPtr v_GetLinStdExp (void) const override
 
void v_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
int v_NumBndryCoeffs () const override
 
int v_NumDGBndryCoeffs () const override
 
void v_ComputeTraceNormal (const int vertex) override
 
void v_ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int mode_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType) override
 Unpack data from input file assuming it comes from. More...
 
const Array< OneD, const NekDouble > & v_GetPhysNormals () override
 
void v_LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
void v_HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
DNekMatSharedPtr v_GenMatrix (const StdRegions::StdMatrixKey &mkey) override
 
DNekScalMatSharedPtr CreateMatrix (const MatrixKey &mkey)
 
DNekMatSharedPtr v_CreateStdMatrix (const StdRegions::StdMatrixKey &mkey) override
 
DNekScalMatSharedPtr v_GetLocMatrix (const MatrixKey &mkey) override
 
void v_DropLocMatrix (const MatrixKey &mkey) override
 
DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const MatrixKey &mkey) override
 
void v_DropLocStaticCondMatrix (const MatrixKey &mkey) override
 
- Protected Member Functions inherited from Nektar::StdRegions::StdSegExp
NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray) override
 Integrate the physical point list inarray over region and return the value. More...
 
void v_PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray) override
 Evaluate the derivative \( d/d{\xi_1} \) at the physical quadrature points given by inarray and return in outarray. 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=NullNekDouble1DArray, 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 transform from coefficient space given in inarray and evaluate at the physical quadrature points outarray. More...
 
void v_BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_FwdTransBndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Forward transform from physical quadrature space stored in inarray and evaluate the expansion coefficients and store in outarray. More...
 
void v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Inner product of inarray over region with respect to the expansion basis (this)->m_base[0] and return in outarray. More...
 
void v_IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check) override
 Inner product of inarray over region with respect to expansion basis base and return in outarray. More...
 
void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) 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
 
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
 
NekDouble v_PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs, std::array< NekDouble, 6 > &secondOrderDerivs) override
 
void v_LaplacianMatrixOp (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_SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey) override
 
void v_ExponentialFilter (Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff) override
 
void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_FillMode (const int mode, Array< OneD, NekDouble > &outarray) override
 
void v_GetCoords (Array< OneD, NekDouble > &coords_0, Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2) override
 
void v_GetBoundaryMap (Array< OneD, unsigned int > &outarray) override
 
void v_GetInteriorMap (Array< OneD, unsigned int > &outarray) override
 
int v_GetVertexMap (int localVertexId, bool useCoeffPacking=false) override
 
int v_GetNverts () const final
 
int v_GetNtraces () const final
 
int v_GetTraceNcoeffs (const int i) const final
 
int v_GetTraceIntNcoeffs (const int i) const final
 
int v_GetTraceNumPoints (const int i) const final
 
int v_NumBndryCoeffs () const override
 
int v_NumDGBndryCoeffs () const override
 
bool v_IsBoundaryInteriorExpansion () const override
 
int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset) override
 
LibUtilities::ShapeType v_DetShapeType () const override
 Return Shape of region, using ShapeType enum list. i.e. Segment. More...
 
DNekMatSharedPtr v_GenMatrix (const StdMatrixKey &mkey) override
 
DNekMatSharedPtr v_CreateStdMatrix (const StdMatrixKey &mkey) override
 
void v_GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true) override
 
void v_ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_GetTraceCoeffMap (const unsigned int traceid, Array< OneD, unsigned int > &maparray) override
 Get the map of the coefficient location to teh local trace coefficients. More...
 
void v_GetElmtTraceToTraceMap (const unsigned int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient, int P, int Q) override
 
void v_GetTraceToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient, int P, int Q) override
 
NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &coords, 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
 
NekDouble v_PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs, std::array< NekDouble, 6 > &secondOrderDerivs) 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)
 
- Protected Member Functions inherited from Nektar::LocalRegions::Expansion1D
DNekMatSharedPtr v_GenMatrix (const StdRegions::StdMatrixKey &mkey) override
 
void v_AddRobinMassMatrix (const int vert, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat) override
 
void v_AddRobinTraceContribution (const int vert, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs) override
 
NekDouble v_VectorFlux (const Array< OneD, Array< OneD, NekDouble > > &vec) override
 
void v_NormalTraceDerivFactors (Array< OneD, Array< OneD, NekDouble > > &factors, Array< OneD, Array< OneD, NekDouble > > &d0factors, Array< OneD, Array< OneD, NekDouble > > &d1factors) override
 : This method gets all of the factors which are required as part of the Gradient Jump Penalty stabilisation and involves the product of the normal and geometric factors along the element trace. More...
 
void v_ReOrientTracePhysMap (const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1) override
 
void v_TraceNormLen (const int traceid, NekDouble &h, NekDouble &p) override
 
- Protected Member Functions inherited from Nektar::LocalRegions::Expansion
void ComputeLaplacianMetric ()
 
void ComputeQuadratureMetric ()
 
void ComputeGmatcdotMF (const Array< TwoD, const NekDouble > &df, const Array< OneD, const NekDouble > &direction, Array< OneD, Array< OneD, NekDouble > > &dfdir)
 
Array< OneD, NekDoubleGetMF (const int dir, const int shapedim, const StdRegions::VarCoeffMap &varcoeffs)
 
Array< OneD, NekDoubleGetMFDiv (const int dir, const StdRegions::VarCoeffMap &varcoeffs)
 
Array< OneD, NekDoubleGetMFMag (const int dir, const StdRegions::VarCoeffMap &varcoeffs)
 
void v_MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
virtual void v_DivideByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_ComputeLaplacianMetric ()
 
int v_GetCoordim () const override
 
void v_GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
 
virtual DNekScalMatSharedPtr v_GetLocMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual void v_DropLocMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual DNekMatSharedPtr v_BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
 
virtual DNekMatSharedPtr v_BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd)
 
virtual void v_ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType)
 
virtual void v_AddEdgeNormBoundaryInt (const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddEdgeNormBoundaryInt (const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddFaceNormBoundaryInt (const int face, const std::shared_ptr< Expansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
virtual void v_DGDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &coeffs, Array< OneD, NekDouble > &outarray)
 
virtual NekDouble v_VectorFlux (const Array< OneD, Array< OneD, NekDouble > > &vec)
 
virtual void v_NormalTraceDerivFactors (Array< OneD, Array< OneD, NekDouble > > &factors, Array< OneD, Array< OneD, NekDouble > > &d0factors, Array< OneD, Array< OneD, NekDouble > > &d1factors)
 
virtual void v_AlignVectorToCollapsedDir (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual StdRegions::Orientation v_GetTraceOrient (int trace)
 
void v_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
virtual void v_GetTraceQFactors (const int trace, Array< OneD, NekDouble > &outarray)
 
virtual void v_GetTracePhysVals (const int trace, const StdRegions::StdExpansionSharedPtr &TraceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient)
 
virtual void v_GetTracePhysMap (const int edge, Array< OneD, int > &outarray)
 
virtual void v_ReOrientTracePhysMap (const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1=-1)
 
virtual void v_ComputeTraceNormal (const int id)
 
virtual const Array< OneD, const NekDouble > & v_GetPhysNormals ()
 
virtual void v_SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
virtual void v_SetUpPhysNormals (const int id)
 
virtual void v_AddRobinMassMatrix (const int face, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
virtual void v_AddRobinTraceContribution (const int traceid, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs)
 
virtual void v_TraceNormLen (const int traceid, NekDouble &h, NekDouble &p)
 
virtual void v_GenTraceExp (const int traceid, ExpansionSharedPtr &exp)
 

Private Member Functions

void ReverseCoeffsAndSign (const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Reverse the coefficients in a boundary interior expansion this routine is of use when we need the segment coefficients corresponding to a expansion in the reverse coordinate direction. More...
 
void MultiplyByElmtInvMass (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 

Private Attributes

LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLessm_matrixManager
 
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLessm_staticCondMatrixManager
 

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
 
- Protected Attributes inherited from Nektar::LocalRegions::Expansion
LibUtilities::NekManager< IndexMapKey, IndexMapValues, IndexMapKey::opLessm_indexMapManager
 
std::map< int, ExpansionWeakPtrm_traceExp
 
SpatialDomains::GeometrySharedPtr m_geom
 
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
 
MetricMap m_metrics
 
std::map< int, NormalVectorm_traceNormals
 
ExpansionWeakPtr m_elementLeft
 
ExpansionWeakPtr m_elementRight
 
int m_elementTraceLeft = -1
 
int m_elementTraceRight = -1
 
std::map< int, Array< OneD, NekDouble > > m_elmtBndNormDirElmtLen
 the element length in each element boundary(Vertex, edge or face) normal direction calculated based on the local m_metricinfo times the standard element length (which is 2.0) More...
 

Detailed Description

Defines a Segment local expansion.

Definition at line 49 of file SegExp.h.

Constructor & Destructor Documentation

◆ SegExp() [1/2]

Nektar::LocalRegions::SegExp::SegExp ( const LibUtilities::BasisKey Ba,
const SpatialDomains::Geometry1DSharedPtr geom 
)

Constructor using BasisKey class for quadrature points and order definition.

Parameters
BaBasis key of segment expansion.
geomDescription of geometry.

Definition at line 55 of file SegExp.cpp.

57 : StdExpansion(Ba.GetNumModes(), 1, Ba),
58 StdExpansion1D(Ba.GetNumModes(), Ba), StdRegions::StdSegExp(Ba),
59 Expansion(geom), Expansion1D(geom),
61 std::bind(&SegExp::CreateMatrix, this, std::placeholders::_1),
62 std::string("SegExpMatrix")),
64 this, std::placeholders::_1),
65 std::string("SegExpStaticCondMatrix"))
66{
67}
Expansion1D(SpatialDomains::Geometry1DSharedPtr pGeom)
Definition: Expansion1D.h:58
DNekScalBlkMatSharedPtr CreateStaticCondMatrix(const MatrixKey &mkey)
Definition: Expansion.cpp:272
Expansion(SpatialDomains::GeometrySharedPtr pGeom)
Definition: Expansion.cpp:43
DNekScalMatSharedPtr CreateMatrix(const MatrixKey &mkey)
Definition: SegExp.cpp:1091
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: SegExp.h:235
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: SegExp.h:233
StdExpansion()
Default Constructor.

◆ SegExp() [2/2]

Nektar::LocalRegions::SegExp::SegExp ( const SegExp S)

Copy Constructor.

Parameters
SExisting segment to duplicate.

Definition at line 73 of file SegExp.cpp.

74 : StdExpansion(S), StdExpansion1D(S), StdRegions::StdSegExp(S),
75 Expansion(S), Expansion1D(S), m_matrixManager(S.m_matrixManager),
76 m_staticCondMatrixManager(S.m_staticCondMatrixManager)
77{
78}

◆ ~SegExp()

Nektar::LocalRegions::SegExp::~SegExp ( )
overridedefault

Member Function Documentation

◆ CreateMatrix()

DNekScalMatSharedPtr Nektar::LocalRegions::SegExp::CreateMatrix ( const MatrixKey mkey)
protected

Definition at line 1091 of file SegExp.cpp.

1092{
1093 DNekScalMatSharedPtr returnval;
1094 NekDouble fac;
1096
1098 "Geometric information is not set up");
1099
1100 switch (mkey.GetMatrixType())
1101 {
1102 case StdRegions::eMass:
1103 {
1104 if ((m_metricinfo->GetGtype() == SpatialDomains::eDeformed) ||
1105 (mkey.GetNVarCoeff()))
1106 {
1107 fac = 1.0;
1108 goto UseLocRegionsMatrix;
1109 }
1110 else
1111 {
1112 fac = (m_metricinfo->GetJac(ptsKeys))[0];
1113 goto UseStdRegionsMatrix;
1114 }
1115 }
1116 break;
1118 {
1119 if ((m_metricinfo->GetGtype() == SpatialDomains::eDeformed) ||
1120 (mkey.GetNVarCoeff()))
1121 {
1122 NekDouble one = 1.0;
1123 StdRegions::StdMatrixKey masskey(StdRegions::eMass,
1124 DetShapeType(), *this);
1125 DNekMatSharedPtr mat = GenMatrix(masskey);
1126 mat->Invert();
1127
1128 returnval =
1130 }
1131 else
1132 {
1133 fac = 1.0 / (m_metricinfo->GetJac(ptsKeys))[0];
1134 goto UseStdRegionsMatrix;
1135 }
1136 }
1137 break;
1141 {
1142 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed ||
1143 mkey.GetNVarCoeff())
1144 {
1145 fac = 1.0;
1146 goto UseLocRegionsMatrix;
1147 }
1148 else
1149 {
1150 int dir = 0;
1151 switch (mkey.GetMatrixType())
1152 {
1154 dir = 0;
1155 break;
1157 ASSERTL1(m_geom->GetCoordim() >= 2,
1158 "Cannot call eWeakDeriv2 in a "
1159 "coordinate system which is not at "
1160 "least two-dimensional");
1161 dir = 1;
1162 break;
1164 ASSERTL1(m_geom->GetCoordim() == 3,
1165 "Cannot call eWeakDeriv2 in a "
1166 "coordinate system which is not "
1167 "three-dimensional");
1168 dir = 2;
1169 break;
1170 default:
1171 break;
1172 }
1173
1174 MatrixKey deriv0key(StdRegions::eWeakDeriv0,
1175 mkey.GetShapeType(), *this);
1176
1177 DNekMatSharedPtr WeakDerivStd = GetStdMatrix(deriv0key);
1178 fac = m_metricinfo->GetDerivFactors(ptsKeys)[dir][0] *
1179 m_metricinfo->GetJac(ptsKeys)[0];
1180
1182 fac, WeakDerivStd);
1183 }
1184 }
1185 break;
1187 {
1188 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1189 {
1190 fac = 1.0;
1191 goto UseLocRegionsMatrix;
1192 }
1193 else
1194 {
1195 int coordim = m_geom->GetCoordim();
1196 fac = 0.0;
1197 for (int i = 0; i < coordim; ++i)
1198 {
1199 fac += m_metricinfo->GetDerivFactors(ptsKeys)[i][0] *
1200 m_metricinfo->GetDerivFactors(ptsKeys)[i][0];
1201 }
1202 fac *= m_metricinfo->GetJac(ptsKeys)[0];
1203 goto UseStdRegionsMatrix;
1204 }
1205 }
1206 break;
1208 {
1209 NekDouble factor = mkey.GetConstFactor(StdRegions::eFactorLambda);
1210 MatrixKey masskey(StdRegions::eMass, mkey.GetShapeType(), *this);
1211 DNekScalMat &MassMat = *(this->m_matrixManager[masskey]);
1212 MatrixKey lapkey(StdRegions::eLaplacian, mkey.GetShapeType(), *this,
1213 mkey.GetConstFactors(), mkey.GetVarCoeffs());
1214 DNekScalMat &LapMat = *(this->m_matrixManager[lapkey]);
1215
1216 int rows = LapMat.GetRows();
1217 int cols = LapMat.GetColumns();
1218
1219 DNekMatSharedPtr helm =
1221
1222 NekDouble one = 1.0;
1223 (*helm) = LapMat + factor * MassMat;
1224
1225 returnval =
1227 }
1228 break;
1233 {
1234 NekDouble one = 1.0;
1235
1236 DNekMatSharedPtr mat = GenMatrix(mkey);
1238 }
1239 break;
1241 {
1242 NekDouble one = 1.0;
1243
1244 // StdRegions::StdMatrixKey
1245 // hkey(StdRegions::eHybridDGHelmholtz,
1246 // DetShapeType(),*this,
1247 // mkey.GetConstant(0),
1248 // mkey.GetConstant(1));
1250 *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
1251 DNekMatSharedPtr mat = GenMatrix(hkey);
1252
1253 mat->Invert();
1255 }
1256 break;
1258 {
1259 DNekMatSharedPtr m_Ix;
1260 Array<OneD, NekDouble> coords(1, 0.0);
1261 StdRegions::ConstFactorMap factors = mkey.GetConstFactors();
1262 int vertex = (int)factors[StdRegions::eFactorGaussVertex];
1263
1264 coords[0] = (vertex == 0) ? -1.0 : 1.0;
1265
1266 m_Ix = m_base[0]->GetI(coords);
1267 returnval =
1269 }
1270 break;
1271
1272 UseLocRegionsMatrix:
1273 {
1274 DNekMatSharedPtr mat = GenMatrix(mkey);
1276 }
1277 break;
1278 UseStdRegionsMatrix:
1279 {
1280 DNekMatSharedPtr mat = GetStdMatrix(mkey);
1282 }
1283 break;
1284 default:
1285 NEKERROR(ErrorUtil::efatal, "Matrix creation not defined");
1286 break;
1287 }
1288
1289 return returnval;
1290}
#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
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:265
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:273
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:274
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:603
const LibUtilities::PointsKeyVector GetPointsKeys() const
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:367
DNekMatSharedPtr GenMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:844
Array< OneD, LibUtilities::BasisSharedPtr > m_base
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:231
@ eNoGeomType
No type defined.
@ eDeformed
Geometry is curved or has non-constant factors.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:402
StdRegions::ConstFactorMap factors
NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag > DNekScalMat
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL1, ASSERTL2, Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::eFactorGaussVertex, Nektar::StdRegions::eFactorLambda, Nektar::ErrorUtil::efatal, Nektar::StdRegions::eHelmholtz, Nektar::StdRegions::eHybridDGHelmBndLam, Nektar::StdRegions::eHybridDGHelmholtz, Nektar::StdRegions::eHybridDGLamToQ0, Nektar::StdRegions::eHybridDGLamToU, Nektar::StdRegions::eInterpGauss, Nektar::StdRegions::eInvHybridDGHelmholtz, Nektar::StdRegions::eInvMass, Nektar::StdRegions::eLaplacian, Nektar::StdRegions::eMass, Nektar::SpatialDomains::eNoGeomType, Nektar::StdRegions::eWeakDeriv0, Nektar::StdRegions::eWeakDeriv1, Nektar::StdRegions::eWeakDeriv2, Nektar::VarcoeffHashingTest::factors, Nektar::StdRegions::StdExpansion::GenMatrix(), Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdMatrixKey::GetConstFactors(), Nektar::StdRegions::StdMatrixKey::GetMatrixType(), Nektar::StdRegions::StdMatrixKey::GetNVarCoeff(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdMatrixKey::GetShapeType(), Nektar::StdRegions::StdExpansion::GetStdMatrix(), Nektar::StdRegions::StdMatrixKey::GetVarCoeffs(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_geom, m_matrixManager, Nektar::LocalRegions::Expansion::m_metricinfo, and NEKERROR.

◆ MultiplyByElmtInvMass()

void Nektar::LocalRegions::SegExp::MultiplyByElmtInvMass ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
private
Todo:
Same method exists in ExpList and everyone references ExpList::MultiplyByElmtInvMass. Remove this one?

Definition at line 1363 of file SegExp.cpp.

1365{
1366 // get Mass matrix inverse
1367 MatrixKey masskey(StdRegions::eInvMass, DetShapeType(), *this);
1368 DNekScalMatSharedPtr matsys = m_matrixManager[masskey];
1369
1370 NekVector<NekDouble> in(m_ncoeffs, inarray, eCopy);
1371 NekVector<NekDouble> out(m_ncoeffs, outarray, eWrapper);
1372
1373 out = (*matsys) * in;
1374}

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

◆ ReverseCoeffsAndSign()

void Nektar::LocalRegions::SegExp::ReverseCoeffsAndSign ( const Array< OneD, NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
private

Reverse the coefficients in a boundary interior expansion this routine is of use when we need the segment coefficients corresponding to a expansion in the reverse coordinate direction.

Definition at line 1322 of file SegExp.cpp.

1324{
1325
1326 int m;
1327 NekDouble sgn = 1;
1328
1329 ASSERTL1(&inarray[0] != &outarray[0],
1330 "inarray and outarray can not be the same");
1331 switch (GetBasisType(0))
1332 {
1334 // Swap vertices
1335 outarray[0] = inarray[1];
1336 outarray[1] = inarray[0];
1337 // negate odd modes
1338 for (m = 2; m < m_ncoeffs; ++m)
1339 {
1340 outarray[m] = sgn * inarray[m];
1341 sgn = -sgn;
1342 }
1343 break;
1346 for (m = 0; m < m_ncoeffs; ++m)
1347 {
1348 outarray[m_ncoeffs - 1 - m] = inarray[m];
1349 }
1350 break;
1351 default:
1352 ASSERTL0(false, "This basis is not allowed in this method");
1353 break;
1354 }
1355}
#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
@ 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::LibUtilities::eGauss_Lagrange, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), and Nektar::StdRegions::StdExpansion::m_ncoeffs.

Referenced by v_SetCoeffsToOrientation().

◆ v_ComputeTraceNormal()

void Nektar::LocalRegions::SegExp::v_ComputeTraceNormal ( const int  vertex)
overrideprotectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 805 of file SegExp.cpp.

806{
807 int i;
808 const SpatialDomains::GeomFactorsSharedPtr &geomFactors =
809 GetGeom()->GetMetricInfo();
810 SpatialDomains::GeomType type = geomFactors->GetGtype();
811 const Array<TwoD, const NekDouble> &gmat =
812 geomFactors->GetDerivFactors(GetPointsKeys());
813 int nqe = 1;
814 int vCoordDim = GetCoordim();
815
816 m_traceNormals[vertex] = Array<OneD, Array<OneD, NekDouble>>(vCoordDim);
817 Array<OneD, Array<OneD, NekDouble>> &normal = m_traceNormals[vertex];
818 for (i = 0; i < vCoordDim; ++i)
819 {
820 normal[i] = Array<OneD, NekDouble>(nqe);
821 }
822
823 size_t nqb = nqe;
824 size_t nbnd = vertex;
825 m_elmtBndNormDirElmtLen[nbnd] = Array<OneD, NekDouble>{nqb, 0.0};
826 Array<OneD, NekDouble> &length = m_elmtBndNormDirElmtLen[nbnd];
827
828 // Regular geometry case
829 if ((type == SpatialDomains::eRegular) ||
831 {
832 NekDouble vert;
833 // Set up normals
834 switch (vertex)
835 {
836 case 0:
837 for (i = 0; i < vCoordDim; ++i)
838 {
839 Vmath::Fill(nqe, -gmat[i][0], normal[i], 1);
840 }
841 break;
842 case 1:
843 for (i = 0; i < vCoordDim; ++i)
844 {
845 Vmath::Fill(nqe, gmat[i][0], normal[i], 1);
846 }
847 break;
848 default:
849 ASSERTL0(false, "point is out of range (point < 2)");
850 }
851
852 // normalise
853 vert = 0.0;
854 for (i = 0; i < vCoordDim; ++i)
855 {
856 vert += normal[i][0] * normal[i][0];
857 }
858 vert = 1.0 / sqrt(vert);
859
860 Vmath::Fill(nqb, vert, length, 1);
861
862 for (i = 0; i < vCoordDim; ++i)
863 {
864 Vmath::Smul(nqe, vert, normal[i], 1, normal[i], 1);
865 }
866 }
867}
std::map< int, NormalVector > m_traceNormals
Definition: Expansion.h:276
std::map< int, Array< OneD, NekDouble > > m_elmtBndNormDirElmtLen
the element length in each element boundary(Vertex, edge or face) normal direction calculated based o...
Definition: Expansion.h:286
SpatialDomains::GeometrySharedPtr GetGeom() const
Definition: Expansion.cpp:167
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition: GeomFactors.h:60
GeomType
Indicates the type of element geometry.
@ eRegular
Geometry is straight-sided with constant geometric factors.
@ eMovingRegular
Currently unused.
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
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.hpp:54
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References ASSERTL0, Nektar::SpatialDomains::eMovingRegular, Nektar::SpatialDomains::eRegular, Vmath::Fill(), Nektar::StdRegions::StdExpansion::GetCoordim(), Nektar::LocalRegions::Expansion::GetGeom(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::LocalRegions::Expansion::m_elmtBndNormDirElmtLen, Nektar::LocalRegions::Expansion::m_traceNormals, Vmath::Smul(), and tinysimd::sqrt().

◆ v_CreateStdMatrix()

DNekMatSharedPtr Nektar::LocalRegions::SegExp::v_CreateStdMatrix ( const StdRegions::StdMatrixKey mkey)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 1082 of file SegExp.cpp.

1083{
1084 LibUtilities::BasisKey bkey = m_base[0]->GetBasisKey();
1087
1088 return tmp->GetStdMatrix(mkey);
1089}
std::shared_ptr< StdSegExp > StdSegExpSharedPtr
Definition: StdSegExp.h:222

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and Nektar::StdRegions::StdExpansion::m_base.

◆ v_DropLocMatrix()

void Nektar::LocalRegions::SegExp::v_DropLocMatrix ( const MatrixKey mkey)
overrideprotectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1077 of file SegExp.cpp.

1078{
1079 m_matrixManager.DeleteObject(mkey);
1080}

References m_matrixManager.

◆ v_DropLocStaticCondMatrix()

void Nektar::LocalRegions::SegExp::v_DropLocStaticCondMatrix ( const MatrixKey mkey)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1067 of file SegExp.cpp.

1068{
1069 m_staticCondMatrixManager.DeleteObject(mkey);
1070}

References m_staticCondMatrixManager.

◆ v_ExtractDataToCoeffs()

void Nektar::LocalRegions::SegExp::v_ExtractDataToCoeffs ( const NekDouble data,
const std::vector< unsigned int > &  nummodes,
const int  mode_offset,
NekDouble coeffs,
std::vector< LibUtilities::BasisType > &  fromType 
)
overrideprotectedvirtual

Unpack data from input file assuming it comes from.

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 763 of file SegExp.cpp.

767{
768 switch (m_base[0]->GetBasisType())
769 {
771 {
772 int fillorder = min((int)nummodes[mode_offset], m_ncoeffs);
773
774 Vmath::Zero(m_ncoeffs, coeffs, 1);
775 Vmath::Vcopy(fillorder, &data[0], 1, &coeffs[0], 1);
776 }
777 break;
779 {
780 // Assume that input is also Gll_Lagrange
781 // but no way to check;
782 LibUtilities::PointsKey f0(nummodes[mode_offset],
784 LibUtilities::PointsKey t0(m_base[0]->GetNumModes(),
786 LibUtilities::Interp1D(f0, data, t0, coeffs);
787 }
788 break;
790 {
791 // Assume that input is also Gauss_Lagrange
792 // but no way to check;
793 LibUtilities::PointsKey f0(nummodes[mode_offset],
795 LibUtilities::PointsKey t0(m_base[0]->GetNumModes(),
797 LibUtilities::Interp1D(f0, data, t0, coeffs);
798 }
799 break;
800 default:
801 ASSERTL0(false, "basis is either not set up or not hierarchicial");
802 }
803}
void Interp1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
this function interpolates a 1D function evaluated at the quadrature points of the basis fbasis0 to ...
Definition: Interp.cpp:47
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
@ eGaussGaussLegendre
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:46
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References ASSERTL0, Nektar::LibUtilities::eGauss_Lagrange, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::LibUtilities::Interp1D(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Vmath::Vcopy(), and Vmath::Zero().

◆ v_FwdTrans()

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

Forward transform from physical quadrature space stored in inarray and evaluate the expansion coefficients and store in outarray.

Perform a forward transform using a Galerkin projection by taking the inner product of the physical points and multiplying by the inverse of the mass matrix using the Solve method of the standard matrix container holding the local mass matrix, i.e. \( {\bf \hat{u}} = {\bf M}^{-1} {\bf I} \) where \( {\bf I}[p] = \int^1_{-1} \phi_p(\xi_1) u(\xi_1) d\xi_1 \)

Inputs:

  • inarray: array of physical quadrature points to be transformed

Outputs:

  • outarray: updated array of expansion coefficients.

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 346 of file SegExp.cpp.

348{
349 if (m_base[0]->Collocation())
350 {
351 Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
352 }
353 else
354 {
355 v_IProductWRTBase(inarray, outarray);
356
357 // get Mass matrix inverse
358 MatrixKey masskey(StdRegions::eInvMass, DetShapeType(), *this);
359 DNekScalMatSharedPtr matsys = m_matrixManager[masskey];
360
361 // copy inarray in case inarray == outarray
362 NekVector<NekDouble> in(m_ncoeffs, outarray, eCopy);
363 NekVector<NekDouble> out(m_ncoeffs, outarray, eWrapper);
364
365 out = (*matsys) * in;
366 }
367}
void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Inner product of inarray over region with respect to the expansion basis (this)->_Base[0] and return ...
Definition: SegExp.cpp:472

References Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::eCopy, Nektar::StdRegions::eInvMass, Nektar::eWrapper, Nektar::StdRegions::StdExpansion::m_base, m_matrixManager, Nektar::StdRegions::StdExpansion::m_ncoeffs, v_IProductWRTBase(), and Vmath::Vcopy().

Referenced by v_FwdTransBndConstrained().

◆ v_FwdTransBndConstrained()

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

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 369 of file SegExp.cpp.

372{
373 if (m_base[0]->Collocation())
374 {
375 Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
376 }
377 else
378 {
379 int nInteriorDofs = m_ncoeffs - 2;
380 int offset = 0;
381
382 switch (m_base[0]->GetBasisType())
383 {
385 {
386 offset = 1;
387 }
388 break;
390 {
391 nInteriorDofs = m_ncoeffs;
392 offset = 0;
393 }
394 break;
397 {
398 ASSERTL1(
399 m_base[0]->GetPointsType() ==
401 m_base[0]->GetPointsType() ==
403 "Cannot use FwdTrans_BndConstrained with these points.");
404 offset = 2;
405 }
406 break;
407 default:
408 ASSERTL0(false, "This type of FwdTrans is not defined"
409 "for this expansion type");
410 }
411
412 fill(outarray.get(), outarray.get() + m_ncoeffs, 0.0);
413
415 {
416
417 outarray[GetVertexMap(0)] = inarray[0];
418 outarray[GetVertexMap(1)] = inarray[m_base[0]->GetNumPoints() - 1];
419
420 if (m_ncoeffs > 2)
421 {
422 // ideally, we would like to have tmp0 to be replaced
423 // by outarray (currently MassMatrixOp does not allow
424 // aliasing)
425 Array<OneD, NekDouble> tmp0(m_ncoeffs);
426 Array<OneD, NekDouble> tmp1(m_ncoeffs);
427
428 StdRegions::StdMatrixKey stdmasskey(StdRegions::eMass,
429 DetShapeType(), *this);
430 MassMatrixOp(outarray, tmp0, stdmasskey);
431 v_IProductWRTBase(inarray, tmp1);
432
433 Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
434
435 // get Mass matrix inverse (only of interior DOF)
436 MatrixKey masskey(StdRegions::eMass, DetShapeType(), *this);
437 DNekScalMatSharedPtr matsys =
438 (m_staticCondMatrixManager[masskey])->GetBlock(1, 1);
439
440 Blas::Dgemv('N', nInteriorDofs, nInteriorDofs, matsys->Scale(),
441 &((matsys->GetOwnedMatrix())->GetPtr())[0],
442 nInteriorDofs, tmp1.get() + offset, 1, 0.0,
443 outarray.get() + offset, 1);
444 }
445 }
446 else
447 {
448 SegExp::v_FwdTrans(inarray, outarray);
449 }
450 }
451}
void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Forward transform from physical quadrature space stored in inarray and evaluate the expansion coeffic...
Definition: SegExp.cpp:346
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:752
int GetVertexMap(const int localVertexId, bool useCoeffPacking=false)
Definition: StdExpansion.h:679
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 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
@ ePolyEvenlySpaced
1D Evenly-spaced points using Lagrange polynomial
Definition: PointsType.h:73
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:49
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 ASSERTL0, ASSERTL1, Nektar::StdRegions::StdExpansion::DetShapeType(), Blas::Dgemv(), Nektar::LibUtilities::eGauss_Lagrange, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGLL_Lagrange, Nektar::StdRegions::eMass, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::LibUtilities::ePolyEvenlySpaced, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::StdRegions::StdExpansion::GetVertexMap(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, m_staticCondMatrixManager, Nektar::StdRegions::StdExpansion::MassMatrixOp(), v_FwdTrans(), v_IProductWRTBase(), Vmath::Vcopy(), and Vmath::Vsub().

◆ v_GenMatrix()

DNekMatSharedPtr Nektar::LocalRegions::SegExp::v_GenMatrix ( const StdRegions::StdMatrixKey mkey)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 1292 of file SegExp.cpp.

1293{
1294 DNekMatSharedPtr returnval;
1295
1296 switch (mkey.GetMatrixType())
1297 {
1304 returnval = Expansion1D::v_GenMatrix(mkey);
1305 break;
1306 default:
1307 returnval = StdSegExp::v_GenMatrix(mkey);
1308 break;
1309 }
1310
1311 return returnval;
1312}
DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey) override
Definition: Expansion1D.cpp:43

References Nektar::StdRegions::eHybridDGHelmBndLam, Nektar::StdRegions::eHybridDGHelmholtz, Nektar::StdRegions::eHybridDGLamToQ0, Nektar::StdRegions::eHybridDGLamToQ1, Nektar::StdRegions::eHybridDGLamToQ2, Nektar::StdRegions::eHybridDGLamToU, Nektar::StdRegions::StdMatrixKey::GetMatrixType(), and Nektar::LocalRegions::Expansion1D::v_GenMatrix().

◆ v_GetCoord()

void Nektar::LocalRegions::SegExp::v_GetCoord ( const Array< OneD, const NekDouble > &  Lcoords,
Array< OneD, NekDouble > &  coords 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 626 of file SegExp.cpp.

628{
629 int i;
630
631 ASSERTL1(Lcoords[0] >= -1.0 && Lcoords[0] <= 1.0,
632 "Local coordinates are not in region [-1,1]");
633
634 m_geom->FillGeom();
635 for (i = 0; i < m_geom->GetCoordim(); ++i)
636 {
637 coords[i] = m_geom->GetCoord(i, Lcoords);
638 }
639}

References ASSERTL1, and Nektar::LocalRegions::Expansion::m_geom.

◆ v_GetCoords()

void Nektar::LocalRegions::SegExp::v_GetCoords ( Array< OneD, NekDouble > &  coords_1,
Array< OneD, NekDouble > &  coords_2,
Array< OneD, NekDouble > &  coords_3 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 641 of file SegExp.cpp.

644{
645 Expansion::v_GetCoords(coords_0, coords_1, coords_2);
646}
void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
Definition: Expansion.cpp:530

References Nektar::LocalRegions::Expansion::v_GetCoords().

◆ v_GetLinStdExp()

StdRegions::StdExpansionSharedPtr Nektar::LocalRegions::SegExp::v_GetLinStdExp ( void  ) const
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 737 of file SegExp.cpp.

738{
739 LibUtilities::BasisKey bkey0(m_base[0]->GetBasisType(), 2,
740 m_base[0]->GetPointsKey());
741
743}

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

◆ v_GetLocMatrix()

DNekScalMatSharedPtr Nektar::LocalRegions::SegExp::v_GetLocMatrix ( const MatrixKey mkey)
overrideprotectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1072 of file SegExp.cpp.

1073{
1074 return m_matrixManager[mkey];
1075}

References m_matrixManager.

◆ v_GetLocStaticCondMatrix()

DNekScalBlkMatSharedPtr Nektar::LocalRegions::SegExp::v_GetLocStaticCondMatrix ( const MatrixKey mkey)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1062 of file SegExp.cpp.

1063{
1064 return m_staticCondMatrixManager[mkey];
1065}

References m_staticCondMatrixManager.

◆ v_GetPhysNormals()

const Array< OneD, const NekDouble > & Nektar::LocalRegions::SegExp::v_GetPhysNormals ( void  )
overrideprotectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 745 of file SegExp.cpp.

746{
747 NEKERROR(ErrorUtil::efatal, "Got to SegExp");
749}
static Array< OneD, NekDouble > NullNekDouble1DArray

References Nektar::ErrorUtil::efatal, NEKERROR, and Nektar::NullNekDouble1DArray.

◆ v_GetStdExp()

StdRegions::StdExpansionSharedPtr Nektar::LocalRegions::SegExp::v_GetStdExp ( void  ) const
overrideprotectedvirtual

◆ v_GetTracePhysMap()

void Nektar::LocalRegions::SegExp::v_GetTracePhysMap ( const int  vertex,
Array< OneD, int > &  map 
)
overrideprotectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 697 of file SegExp.cpp.

698{
699 int nquad = m_base[0]->GetNumPoints();
700
701 ASSERTL1(vertex == 0 || vertex == 1, "Vertex value should be 0 or 1");
702
703 map = Array<OneD, int>(1);
704
705 map[0] = vertex == 0 ? 0 : nquad - 1;
706}

References ASSERTL1, and Nektar::StdRegions::StdExpansion::m_base.

◆ v_GetTracePhysVals()

void Nektar::LocalRegions::SegExp::v_GetTracePhysVals ( const int  edge,
const StdRegions::StdExpansionSharedPtr EdgeExp,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
StdRegions::Orientation  orient 
)
overrideprotectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 684 of file SegExp.cpp.

690{
691 NekDouble result;
692 v_GetVertexPhysVals(edge, inarray, result);
693 outarray[0] = result;
694}
void v_GetVertexPhysVals(const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray) override
Definition: SegExp.cpp:649

References v_GetVertexPhysVals().

◆ v_GetVertexPhysVals()

void Nektar::LocalRegions::SegExp::v_GetVertexPhysVals ( const int  vertex,
const Array< OneD, const NekDouble > &  inarray,
NekDouble outarray 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 649 of file SegExp.cpp.

652{
653 int nquad = m_base[0]->GetNumPoints();
654
656 {
657 switch (vertex)
658 {
659 case 0:
660 outarray = inarray[0];
661 break;
662 case 1:
663 outarray = inarray[nquad - 1];
664 break;
665 }
666 }
667 else
668 {
671
672 StdRegions::StdMatrixKey key(StdRegions::eInterpGauss, DetShapeType(),
673 *this, factors);
674
675 DNekScalMatSharedPtr mat_gauss = m_matrixManager[key];
676
677 outarray =
678 Blas::Ddot(nquad, mat_gauss->GetOwnedMatrix()->GetPtr().get(), 1,
679 &inarray[0], 1);
680 }
681}
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(), Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::StdRegions::eFactorGaussVertex, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::StdRegions::eInterpGauss, Nektar::VarcoeffHashingTest::factors, Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::StdRegions::StdExpansion::m_base, and m_matrixManager.

Referenced by v_GetTracePhysVals().

◆ v_HelmholtzMatrixOp()

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

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 963 of file SegExp.cpp.

966{
967 int nquad = m_base[0]->GetNumPoints();
968 const Array<TwoD, const NekDouble> &gmat =
969 m_metricinfo->GetDerivFactors(GetPointsKeys());
970 const NekDouble lambda = mkey.GetConstFactor(StdRegions::eFactorLambda);
971
972 Array<OneD, NekDouble> physValues(nquad);
973 Array<OneD, NekDouble> dPhysValuesdx(nquad);
974 Array<OneD, NekDouble> wsp(m_ncoeffs);
975
976 BwdTrans(inarray, physValues);
977
978 // mass matrix operation
979 v_IProductWRTBase((m_base[0]->GetBdata()), physValues, wsp, 1);
980
981 // Laplacian matrix operation
982 switch (m_geom->GetCoordim())
983 {
984 case 1:
985 {
986 PhysDeriv(physValues, dPhysValuesdx);
987
988 // multiply with the proper geometric factors
989 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
990 {
991 Vmath::Vmul(nquad, &gmat[0][0], 1, dPhysValuesdx.get(), 1,
992 dPhysValuesdx.get(), 1);
993 }
994 else
995 {
996 Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.get(), 1);
997 }
998 }
999 break;
1000 case 2:
1001 {
1002 Array<OneD, NekDouble> dPhysValuesdy(nquad);
1003
1004 PhysDeriv(physValues, dPhysValuesdx, dPhysValuesdy);
1005
1006 // multiply with the proper geometric factors
1007 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1008 {
1009 Vmath::Vmul(nquad, &gmat[0][0], 1, dPhysValuesdx.get(), 1,
1010 dPhysValuesdx.get(), 1);
1011 Vmath::Vvtvp(nquad, &gmat[1][0], 1, dPhysValuesdy.get(), 1,
1012 dPhysValuesdx.get(), 1, dPhysValuesdx.get(), 1);
1013 }
1014 else
1015 {
1016 Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.get(), 1);
1017 Blas::Daxpy(nquad, gmat[1][0], dPhysValuesdy.get(), 1,
1018 dPhysValuesdx.get(), 1);
1019 }
1020 }
1021 break;
1022 case 3:
1023 {
1024 Array<OneD, NekDouble> dPhysValuesdy(nquad);
1025 Array<OneD, NekDouble> dPhysValuesdz(nquad);
1026
1027 PhysDeriv(physValues, dPhysValuesdx, dPhysValuesdy, dPhysValuesdz);
1028
1029 // multiply with the proper geometric factors
1030 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1031 {
1032 Vmath::Vmul(nquad, &gmat[0][0], 1, dPhysValuesdx.get(), 1,
1033 dPhysValuesdx.get(), 1);
1034 Vmath::Vvtvp(nquad, &gmat[1][0], 1, dPhysValuesdy.get(), 1,
1035 dPhysValuesdx.get(), 1, dPhysValuesdx.get(), 1);
1036 Vmath::Vvtvp(nquad, &gmat[2][0], 1, dPhysValuesdz.get(), 1,
1037 dPhysValuesdx.get(), 1, dPhysValuesdx.get(), 1);
1038 }
1039 else
1040 {
1041 Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.get(), 1);
1042 Blas::Daxpy(nquad, gmat[1][0], dPhysValuesdy.get(), 1,
1043 dPhysValuesdx.get(), 1);
1044 Blas::Daxpy(nquad, gmat[2][0], dPhysValuesdz.get(), 1,
1045 dPhysValuesdx.get(), 1);
1046 }
1047 }
1048 break;
1049 default:
1050 ASSERTL0(false, "Wrong number of dimensions");
1051 break;
1052 }
1053
1054 v_IProductWRTBase(m_base[0]->GetDbdata(), dPhysValuesdx, outarray, 1);
1055 Blas::Daxpy(m_ncoeffs, lambda, wsp.get(), 1, outarray.get(), 1);
1056}
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Backward transformation from coefficient space to physical space.
Definition: StdExpansion.h:424
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
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition: Blas.hpp:149
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
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 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

References ASSERTL0, Nektar::StdRegions::StdExpansion::BwdTrans(), Blas::Daxpy(), Blas::Dscal(), Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_geom, Nektar::LocalRegions::Expansion::m_metricinfo, Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::StdRegions::StdExpansion::PhysDeriv(), v_IProductWRTBase(), Vmath::Vmul(), and Vmath::Vvtvp().

◆ v_Integral()

NekDouble Nektar::LocalRegions::SegExp::v_Integral ( const Array< OneD, const NekDouble > &  inarray)
overrideprotectedvirtual

Integrate the physical point list inarray over region and return the value.

Inputs:

  • inarray: definition of function to be returned at quadrature point of expansion.

Outputs:

  • returns \(\int^1_{-1} u(\xi_1)d \xi_1 \) where \(inarray[i] = u(\xi_{1i}) \)

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 98 of file SegExp.cpp.

99{
100 int nquad0 = m_base[0]->GetNumPoints();
101 Array<OneD, const NekDouble> jac = m_metricinfo->GetJac(GetPointsKeys());
102 NekDouble ival;
103 Array<OneD, NekDouble> tmp(nquad0);
104
105 // multiply inarray with Jacobian
106 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
107 {
108 Vmath::Vmul(nquad0, jac, 1, inarray, 1, tmp, 1);
109 }
110 else
111 {
112 Vmath::Smul(nquad0, jac[0], inarray, 1, tmp, 1);
113 }
114
115 // call StdSegExp version;
116 ival = StdSegExp::v_Integral(tmp);
117 // ival = StdSegExp::Integral(tmp);
118 return ival;
119}

References Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_metricinfo, Vmath::Smul(), and Vmath::Vmul().

◆ v_IProductWRTBase() [1/2]

void Nektar::LocalRegions::SegExp::v_IProductWRTBase ( const Array< OneD, const NekDouble > &  base,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
int  coll_check 
)
overrideprotectedvirtual

Inner product of inarray over region with respect to expansion basis base and return in outarray.

Calculate \( I[p] = \int^{1}_{-1} \phi_p(\xi_1) u(\xi_1) d\xi_1 = \sum_{i=0}^{nq-1} \phi_p(\xi_{1i}) u(\xi_{1i}) w_i \) where \( outarray[p] = I[p], inarray[i] = u(\xi_{1i}), base[p*nq+i] = \phi_p(\xi_{1i}) \).

Inputs:

  • base: an array definiing the local basis for the inner product usually passed from Basis->get_bdata() or Basis->get_Dbdata()
  • inarray: physical point array of function to be integrated \( u(\xi_1) \)
  • coll_check: Flag to identify when a Basis->collocation() call should be performed to see if this is a GLL_Lagrange basis with a collocation property. (should be set to 0 if taking the inner product with respect to the derivative of basis)

Output:

  • outarray: array of coefficients representing the inner product of function with ever mode in the exapnsion

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 505 of file SegExp.cpp.

508{
509 int nquad0 = m_base[0]->GetNumPoints();
510 Array<OneD, const NekDouble> jac = m_metricinfo->GetJac(GetPointsKeys());
511 Array<OneD, NekDouble> tmp(nquad0);
512
513 // multiply inarray with Jacobian
514 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
515 {
516 Vmath::Vmul(nquad0, jac, 1, inarray, 1, tmp, 1);
517 }
518 else
519 {
520 Vmath::Smul(nquad0, jac[0], inarray, 1, tmp, 1);
521 }
522 StdSegExp::v_IProductWRTBase(base, tmp, outarray, coll_check);
523}

References Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_metricinfo, Vmath::Smul(), and Vmath::Vmul().

◆ v_IProductWRTBase() [2/2]

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

Inner product of inarray over region with respect to the expansion basis (this)->_Base[0] and return in outarray.

Wrapper call to SegExp::IProduct_WRT_B

Input:

  • inarray: array of function evaluated at the physical collocation points

Output:

  • outarray: array of inner product with respect to each basis over region

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 472 of file SegExp.cpp.

474{
475 v_IProductWRTBase(m_base[0]->GetBdata(), inarray, outarray, 1);
476}

References Nektar::StdRegions::StdExpansion::m_base, and v_IProductWRTBase().

Referenced by v_FwdTrans(), v_FwdTransBndConstrained(), v_HelmholtzMatrixOp(), v_IProductWRTBase(), v_IProductWRTDerivBase(), v_LaplacianMatrixOp(), and v_NormVectorIProductWRTBase().

◆ v_IProductWRTDerivBase()

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

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 525 of file SegExp.cpp.

528{
529 ASSERTL1(dir < 3, "input dir is out of range");
530 ASSERTL1((dir == 2) ? m_geom->GetCoordim() == 3 : true,
531 "input dir is out of range");
532
533 int nquad = m_base[0]->GetNumPoints();
534 const Array<TwoD, const NekDouble> &gmat =
535 m_metricinfo->GetDerivFactors(GetPointsKeys());
536
537 Array<OneD, NekDouble> tmp1(nquad);
538
539 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
540 {
541 Vmath::Vmul(nquad, gmat[dir], 1, inarray, 1, tmp1, 1);
542 }
543 else
544 {
545 Vmath::Smul(nquad, gmat[dir][0], inarray, 1, tmp1, 1);
546 }
547
548 v_IProductWRTBase(m_base[0]->GetDbdata(), tmp1, outarray, 1);
549}

References ASSERTL1, Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_geom, Nektar::LocalRegions::Expansion::m_metricinfo, Vmath::Smul(), v_IProductWRTBase(), and Vmath::Vmul().

◆ v_LaplacianMatrixOp()

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

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 873 of file SegExp.cpp.

877{
878 int nquad = m_base[0]->GetNumPoints();
879 const Array<TwoD, const NekDouble> &gmat =
880 m_metricinfo->GetDerivFactors(GetPointsKeys());
881
882 Array<OneD, NekDouble> physValues(nquad);
883 Array<OneD, NekDouble> dPhysValuesdx(nquad);
884
885 BwdTrans(inarray, physValues);
886
887 // Laplacian matrix operation
888 switch (m_geom->GetCoordim())
889 {
890 case 1:
891 {
892 PhysDeriv(physValues, dPhysValuesdx);
893
894 // multiply with the proper geometric factors
895 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
896 {
897 Vmath::Vmul(nquad, &gmat[0][0], 1, dPhysValuesdx.get(), 1,
898 dPhysValuesdx.get(), 1);
899 }
900 else
901 {
902 Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.get(), 1);
903 }
904 }
905 break;
906 case 2:
907 {
908 Array<OneD, NekDouble> dPhysValuesdy(nquad);
909
910 PhysDeriv(physValues, dPhysValuesdx, dPhysValuesdy);
911
912 // multiply with the proper geometric factors
913 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
914 {
915 Vmath::Vmul(nquad, &gmat[0][0], 1, dPhysValuesdx.get(), 1,
916 dPhysValuesdx.get(), 1);
917 Vmath::Vvtvp(nquad, &gmat[1][0], 1, dPhysValuesdy.get(), 1,
918 dPhysValuesdx.get(), 1, dPhysValuesdx.get(), 1);
919 }
920 else
921 {
922 Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.get(), 1);
923 Blas::Daxpy(nquad, gmat[1][0], dPhysValuesdy.get(), 1,
924 dPhysValuesdx.get(), 1);
925 }
926 }
927 break;
928 case 3:
929 {
930 Array<OneD, NekDouble> dPhysValuesdy(nquad);
931 Array<OneD, NekDouble> dPhysValuesdz(nquad);
932
933 PhysDeriv(physValues, dPhysValuesdx, dPhysValuesdy, dPhysValuesdz);
934
935 // multiply with the proper geometric factors
936 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
937 {
938 Vmath::Vmul(nquad, &gmat[0][0], 1, dPhysValuesdx.get(), 1,
939 dPhysValuesdx.get(), 1);
940 Vmath::Vvtvp(nquad, &gmat[1][0], 1, dPhysValuesdy.get(), 1,
941 dPhysValuesdx.get(), 1, dPhysValuesdx.get(), 1);
942 Vmath::Vvtvp(nquad, &gmat[2][0], 1, dPhysValuesdz.get(), 1,
943 dPhysValuesdx.get(), 1, dPhysValuesdx.get(), 1);
944 }
945 else
946 {
947 Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.get(), 1);
948 Blas::Daxpy(nquad, gmat[1][0], dPhysValuesdy.get(), 1,
949 dPhysValuesdx.get(), 1);
950 Blas::Daxpy(nquad, gmat[2][0], dPhysValuesdz.get(), 1,
951 dPhysValuesdx.get(), 1);
952 }
953 }
954 break;
955 default:
956 ASSERTL0(false, "Wrong number of dimensions");
957 break;
958 }
959
960 v_IProductWRTBase(m_base[0]->GetDbdata(), dPhysValuesdx, outarray, 1);
961}

References ASSERTL0, Nektar::StdRegions::StdExpansion::BwdTrans(), Blas::Daxpy(), Blas::Dscal(), Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_geom, Nektar::LocalRegions::Expansion::m_metricinfo, Nektar::StdRegions::StdExpansion::PhysDeriv(), v_IProductWRTBase(), Vmath::Vmul(), and Vmath::Vvtvp().

◆ v_NormVectorIProductWRTBase() [1/2]

void Nektar::LocalRegions::SegExp::v_NormVectorIProductWRTBase ( const Array< OneD, const Array< OneD, NekDouble > > &  Fvec,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 569 of file SegExp.cpp.

572{
573 NormVectorIProductWRTBase(Fvec[0], Fvec[1], outarray);
574}
void NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:613

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

◆ v_NormVectorIProductWRTBase() [2/2]

void Nektar::LocalRegions::SegExp::v_NormVectorIProductWRTBase ( const Array< OneD, const NekDouble > &  Fx,
const Array< OneD, const NekDouble > &  Fy,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 551 of file SegExp.cpp.

554{
555 int nq = m_base[0]->GetNumPoints();
556 Array<OneD, NekDouble> Fn(nq);
557
558 // @TODO: This routine no longer makes sense as a normal is not unique to an
559 // edge
560 const Array<OneD, const Array<OneD, NekDouble>> &normals =
561 GetLeftAdjacentElementExp()->GetTraceNormal(
563 Vmath::Vmul(nq, &Fx[0], 1, &normals[0][0], 1, &Fn[0], 1);
564 Vmath::Vvtvp(nq, &Fy[0], 1, &normals[1][0], 1, &Fn[0], 1, &Fn[0], 1);
565
566 v_IProductWRTBase(Fn, outarray);
567}
ExpansionSharedPtr GetLeftAdjacentElementExp() const
Definition: Expansion.h:441
int GetLeftAdjacentElementTrace() const
Definition: Expansion.h:454

References Nektar::LocalRegions::Expansion::GetLeftAdjacentElementExp(), Nektar::LocalRegions::Expansion::GetLeftAdjacentElementTrace(), Nektar::StdRegions::StdExpansion::m_base, v_IProductWRTBase(), Vmath::Vmul(), and Vmath::Vvtvp().

◆ v_NumBndryCoeffs()

int Nektar::LocalRegions::SegExp::v_NumBndryCoeffs ( ) const
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 751 of file SegExp.cpp.

752{
753 return 2;
754}

◆ v_NumDGBndryCoeffs()

int Nektar::LocalRegions::SegExp::v_NumDGBndryCoeffs ( ) const
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 756 of file SegExp.cpp.

757{
758 return 2;
759}

◆ v_PhysDeriv() [1/2]

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

Evaluate the derivative \( d/d{\xi_1} \) at the physical quadrature points given by inarray and return in outarray.

This is a wrapper around StdExpansion1D::Tensor_Deriv

Input:

  • n: number of derivatives to be evaluated where \( n \leq dim\)
  • inarray: array of function evaluated at the quadrature points

Output:

  • outarray: array of the derivatives \( du/d_{\xi_1}|_{\xi_{1i}} d\xi_1/dx, du/d_{\xi_1}|_{\xi_{1i}} d\xi_1/dy, du/d_{\xi_1}|_{\xi_{1i}} d\xi_1/dz, \) depending on value of dim

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 145 of file SegExp.cpp.

149{
150 int nquad0 = m_base[0]->GetNumPoints();
151 Array<TwoD, const NekDouble> gmat =
152 m_metricinfo->GetDerivFactors(GetPointsKeys());
153 Array<OneD, NekDouble> diff(nquad0);
154
155 // StdExpansion1D::PhysTensorDeriv(inarray,diff);
156 PhysTensorDeriv(inarray, diff);
157 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
158 {
159 if (out_d0.size())
160 {
161 Vmath::Vmul(nquad0, &gmat[0][0], 1, &diff[0], 1, &out_d0[0], 1);
162 }
163
164 if (out_d1.size())
165 {
166 Vmath::Vmul(nquad0, &gmat[1][0], 1, &diff[0], 1, &out_d1[0], 1);
167 }
168
169 if (out_d2.size())
170 {
171 Vmath::Vmul(nquad0, &gmat[2][0], 1, &diff[0], 1, &out_d2[0], 1);
172 }
173 }
174 else
175 {
176 if (out_d0.size())
177 {
178 Vmath::Smul(nquad0, gmat[0][0], diff, 1, out_d0, 1);
179 }
180
181 if (out_d1.size())
182 {
183 Vmath::Smul(nquad0, gmat[1][0], diff, 1, out_d1, 1);
184 }
185
186 if (out_d2.size())
187 {
188 Vmath::Smul(nquad0, gmat[2][0], diff, 1, out_d2, 1);
189 }
190 }
191}
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Evaluate the derivative at the physical quadrature points given by inarray and return in outarray.

References Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_metricinfo, Nektar::StdRegions::StdExpansion1D::PhysTensorDeriv(), Vmath::Smul(), and Vmath::Vmul().

Referenced by v_PhysDeriv_n().

◆ v_PhysDeriv() [2/2]

void Nektar::LocalRegions::SegExp::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::StdSegExp.

Definition at line 289 of file SegExp.cpp.

292{
293 switch (dir)
294 {
295 case 0:
296 {
297 PhysDeriv(inarray, outarray, NullNekDouble1DArray,
299 }
300 break;
301 case 1:
302 {
303 PhysDeriv(inarray, NullNekDouble1DArray, outarray,
305 }
306 break;
307 case 2:
308 {
310 outarray);
311 }
312 break;
313 default:
314 {
315 ASSERTL1(false, "input dir is out of range");
316 }
317 break;
318 }
319}

References ASSERTL1, Nektar::NullNekDouble1DArray, and Nektar::StdRegions::StdExpansion::PhysDeriv().

◆ v_PhysDeriv_n()

void Nektar::LocalRegions::SegExp::v_PhysDeriv_n ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_dn 
)
overrideprotectedvirtual

Evaluate the derivative normal to a line: \( d/dn=\frac{spacedim}{||normal||}d/d{\xi} \). The derivative is calculated performing the product \( du/d{s}=\nabla u \cdot normal \).

Parameters
inarrayfunction to derive
out_dnresult of the derivative operation

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 240 of file SegExp.cpp.

242{
243 int nquad0 = m_base[0]->GetNumPoints();
244 Array<TwoD, const NekDouble> gmat =
245 m_metricinfo->GetDerivFactors(GetPointsKeys());
246 int coordim = m_geom->GetCoordim();
247 Array<OneD, NekDouble> out_dn_tmp(nquad0, 0.0);
248 switch (coordim)
249 {
250 case 2:
251
252 Array<OneD, NekDouble> inarray_d0(nquad0);
253 Array<OneD, NekDouble> inarray_d1(nquad0);
254
255 v_PhysDeriv(inarray, inarray_d0, inarray_d1);
256 Array<OneD, Array<OneD, NekDouble>> normals;
257 normals = Array<OneD, Array<OneD, NekDouble>>(coordim);
258 cout << "der_n" << endl;
259 for (int k = 0; k < coordim; ++k)
260 {
261 normals[k] = Array<OneD, NekDouble>(nquad0);
262 }
263 // @TODO: this routine no longer makes sense, since normals are not
264 // unique on
265 // an edge
266 // normals = GetMetricInfo()->GetNormal();
267 for (int i = 0; i < nquad0; i++)
268 {
269 cout << "nx= " << normals[0][i] << " ny=" << normals[1][i]
270 << endl;
271 }
273 "normal vectors do not exist: check if a"
274 "boundary region is defined as I ");
275 // \nabla u \cdot normal
276 Vmath::Vmul(nquad0, normals[0], 1, inarray_d0, 1, out_dn_tmp, 1);
277 Vmath::Vadd(nquad0, out_dn_tmp, 1, out_dn, 1, out_dn, 1);
278 Vmath::Zero(nquad0, out_dn_tmp, 1);
279 Vmath::Vmul(nquad0, normals[1], 1, inarray_d1, 1, out_dn_tmp, 1);
280 Vmath::Vadd(nquad0, out_dn_tmp, 1, out_dn, 1, out_dn, 1);
281
282 for (int i = 0; i < nquad0; i++)
283 {
284 cout << "deps/dx =" << inarray_d0[i]
285 << " deps/dy=" << inarray_d1[i] << endl;
286 }
287 }
288}
void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray) override
Evaluate the derivative at the physical quadrature points given by inarray and return in outarray.
Definition: SegExp.cpp:145
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
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 ASSERTL0, Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_geom, Nektar::LocalRegions::Expansion::m_metricinfo, Nektar::NullNekDoubleArrayOfArray, v_PhysDeriv(), Vmath::Vadd(), Vmath::Vmul(), and Vmath::Zero().

◆ v_PhysDeriv_s()

void Nektar::LocalRegions::SegExp::v_PhysDeriv_s ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_ds 
)
overrideprotectedvirtual

Evaluate the derivative along a line: \( d/ds=\frac{spacedim}{||tangent||}d/d{\xi} \). The derivative is calculated performing the product \( du/d{s}=\nabla u \cdot tangent \).

Parameters
inarrayfunction to derive
out_dsresult of the derivative operation

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 201 of file SegExp.cpp.

203{
204 int nquad0 = m_base[0]->GetNumPoints();
205 int coordim = m_geom->GetCoordim();
206 Array<OneD, NekDouble> diff(nquad0);
207 // this operation is needed if you put out_ds==inarray
208 Vmath::Zero(nquad0, out_ds, 1);
209 switch (coordim)
210 {
211 case 2:
212 // diff= dU/de
213 Array<OneD, NekDouble> diff(nquad0);
214
215 PhysTensorDeriv(inarray, diff);
216
217 // get dS/de= (Jac)^-1
218 Array<OneD, NekDouble> Jac = m_metricinfo->GetJac(GetPointsKeys());
219 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
220 {
221 // calculate the derivative as (dU/de)*(Jac)^-1
222 Vmath::Vdiv(nquad0, diff, 1, Jac, 1, out_ds, 1);
223 }
224 else
225 {
226 NekDouble invJac = 1 / Jac[0];
227 Vmath::Smul(nquad0, invJac, diff, 1, out_ds, 1);
228 }
229 }
230}
void Vdiv(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:126

References Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_geom, Nektar::LocalRegions::Expansion::m_metricinfo, Nektar::StdRegions::StdExpansion1D::PhysTensorDeriv(), Vmath::Smul(), Vmath::Vdiv(), and Vmath::Zero().

◆ v_PhysEvaluate() [1/3]

NekDouble Nektar::LocalRegions::SegExp::v_PhysEvaluate ( const Array< OneD, const NekDouble > &  coord,
const Array< OneD, const NekDouble > &  physvals 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion1D.

Definition at line 593 of file SegExp.cpp.

595{
596 Array<OneD, NekDouble> Lcoord = Array<OneD, NekDouble>(1);
597
598 ASSERTL0(m_geom, "m_geom not defined");
599 m_geom->GetLocCoords(coord, Lcoord);
600
601 return StdExpansion1D::v_PhysEvaluate(Lcoord, physvals);
602}

References ASSERTL0, and Nektar::LocalRegions::Expansion::m_geom.

◆ v_PhysEvaluate() [2/3]

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

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 604 of file SegExp.cpp.

607{
608 Array<OneD, NekDouble> Lcoord(1);
609 ASSERTL0(m_geom, "m_geom not defined");
610 m_geom->GetLocCoords(coord, Lcoord);
611 return StdSegExp::v_PhysEvaluate(Lcoord, inarray, firstOrderDerivs);
612}

References ASSERTL0, and Nektar::LocalRegions::Expansion::m_geom.

◆ v_PhysEvaluate() [3/3]

NekDouble Nektar::LocalRegions::SegExp::v_PhysEvaluate ( const Array< OneD, NekDouble > &  coord,
const Array< OneD, const NekDouble > &  inarray,
std::array< NekDouble, 3 > &  firstOrderDerivs,
std::array< NekDouble, 6 > &  secondOrderDerivs 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 614 of file SegExp.cpp.

618{
619 Array<OneD, NekDouble> Lcoord(1);
620 ASSERTL0(m_geom, "m_geom not defined");
621 m_geom->GetLocCoords(coord, Lcoord);
622 return StdSegExp::v_PhysEvaluate(Lcoord, inarray, firstOrderDerivs,
623 secondOrderDerivs);
624}

References ASSERTL0, and Nektar::LocalRegions::Expansion::m_geom.

◆ v_SetCoeffsToOrientation()

void Nektar::LocalRegions::SegExp::v_SetCoeffsToOrientation ( StdRegions::Orientation  dir,
Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 712 of file SegExp.cpp.

715{
716
717 if (dir == StdRegions::eBackwards)
718 {
719 if (&inarray[0] != &outarray[0])
720 {
721 Array<OneD, NekDouble> intmp(inarray);
722 ReverseCoeffsAndSign(intmp, outarray);
723 }
724 else
725 {
726 ReverseCoeffsAndSign(inarray, outarray);
727 }
728 }
729}
void ReverseCoeffsAndSign(const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Reverse the coefficients in a boundary interior expansion this routine is of use when we need the seg...
Definition: SegExp.cpp:1322

References Nektar::StdRegions::eBackwards, and ReverseCoeffsAndSign().

◆ v_StdPhysEvaluate()

NekDouble Nektar::LocalRegions::SegExp::v_StdPhysEvaluate ( const Array< OneD, const NekDouble > &  Lcoord,
const Array< OneD, const NekDouble > &  physvals 
)
overrideprotectedvirtual

Given the local cartesian coordinate Lcoord evaluate the value of physvals at this point by calling through to the StdExpansion method

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 585 of file SegExp.cpp.

588{
589 // Evaluate point in local (eta) coordinates.
590 return StdExpansion1D::v_PhysEvaluate(Lcoord, physvals);
591}

Member Data Documentation

◆ m_matrixManager

LibUtilities::NekManager<MatrixKey, DNekScalMat, MatrixKey::opLess> Nektar::LocalRegions::SegExp::m_matrixManager
private

◆ m_staticCondMatrixManager

LibUtilities::NekManager<MatrixKey, DNekScalBlkMat, MatrixKey::opLess> Nektar::LocalRegions::SegExp::m_staticCondMatrixManager
private