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

virtual NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray) override
 Integrate the physical point list inarray over region and return the value. More...
 
virtual 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...
 
virtual void v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Calculate the derivative of the physical points in a given direction. More...
 
virtual void v_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...
 
virtual 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...
 
virtual 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...
 
virtual void v_FwdTransBndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
virtual 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...
 
virtual 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...
 
virtual void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray) override
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble >> &Fvec, Array< OneD, NekDouble > &outarray) override
 
virtual NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals) override
 
virtual NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &coord, const Array< OneD, const NekDouble > &physvals) override
 
virtual void v_GetCoord (const Array< OneD, const NekDouble > &Lcoords, Array< OneD, NekDouble > &coords) override
 
virtual void v_GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
 
virtual void v_GetVertexPhysVals (const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray) override
 
virtual void v_AddVertexPhysVals (const int vertex, const NekDouble &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_GetTracePhysVals (const int edge, const StdRegions::StdExpansionSharedPtr &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient) override
 
virtual void v_GetTracePhysMap (const int vertex, Array< OneD, int > &map) override
 
virtual StdRegions::StdExpansionSharedPtr v_GetStdExp (void) const override
 
virtual StdRegions::StdExpansionSharedPtr v_GetLinStdExp (void) const override
 
virtual int v_GetCoordim () override
 
virtual void v_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
virtual int v_GetNumPoints (const int dir) const
 
virtual int v_GetNcoeffs (void) const
 
virtual const LibUtilities::BasisSharedPtrv_GetBasis (int dir) const
 
virtual int v_NumBndryCoeffs () const override
 
virtual int v_NumDGBndryCoeffs () const override
 
virtual void v_ComputeTraceNormal (const int vertex) override
 
virtual SpatialDomains::GeomType v_MetricInfoType ()
 
virtual 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...
 
virtual const Array< OneD, const NekDouble > & v_GetPhysNormals (void) override
 
virtual void v_LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
virtual void v_HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
virtual DNekMatSharedPtr v_GenMatrix (const StdRegions::StdMatrixKey &mkey) override
 
DNekScalMatSharedPtr CreateMatrix (const MatrixKey &mkey)
 
virtual DNekMatSharedPtr v_CreateStdMatrix (const StdRegions::StdMatrixKey &mkey) override
 
virtual DNekScalMatSharedPtr v_GetLocMatrix (const MatrixKey &mkey) override
 
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const MatrixKey &mkey) override
 
void v_DropLocStaticCondMatrix (const MatrixKey &mkey) override
 
- Protected Member Functions inherited from Nektar::StdRegions::StdSegExp
virtual 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)
 
virtual void v_StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Backward transform from coefficient space given in inarray and evaluate at the physical quadrature points outarray. More...
 
virtual void v_BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
virtual void v_IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 
virtual void v_LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi)
 
NekDouble v_PhysEvaluateBasis (const Array< OneD, const NekDouble > &coords, int mode) final
 
virtual void v_SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
virtual void v_ExponentialFilter (Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
 
virtual void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 
virtual void v_GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
virtual void v_GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
virtual int v_GetVertexMap (int localVertexId, bool useCoeffPacking=false)
 
virtual int v_GetNverts () const
 
virtual int v_GetNtraces () const
 
virtual int v_GetTraceNcoeffs (const int i) const
 
virtual int v_GetTraceNumPoints (const int i) const
 
virtual bool v_IsBoundaryInteriorExpansion ()
 
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
virtual LibUtilities::ShapeType v_DetShapeType () const
 Return Shape of region, using ShapeType enum list. i.e. Segment. More...
 
virtual void v_GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true)
 
virtual void v_ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_GetTraceCoeffMap (const unsigned int traceid, Array< OneD, unsigned int > &maparray)
 Get the map of the coefficient location to teh local trace coefficients. More...
 
virtual void v_GetElmtTraceToTraceMap (const unsigned int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient, int P, int Q)
 
virtual void v_GetTraceToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient, int P, int Q)
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion
DNekMatSharedPtr CreateStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr CreateStdStaticCondMatrix (const StdMatrixKey &mkey)
 Create the static condensation of a matrix when using a boundary interior decomposition. More...
 
void BwdTrans_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase_SumFac (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GeneralMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
 
void LaplacianMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp_MatFree (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void HelmholtzMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
 
virtual void v_GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 
template<int DIR>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals)
 This function performs the barycentric interpolation of the polynomial stored in coord at a point physvals using barycentric interpolation weights in direction. More...
 
template<int DIR>
NekDouble BaryEvaluateBasis (const NekDouble &coord, const int &mode)
 
- Protected Member Functions inherited from Nektar::LocalRegions::Expansion1D
virtual void v_AddRobinMassMatrix (const int vert, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
virtual void v_AddRobinEdgeContribution (const int vert, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs)
 
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)
 : 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...
 
virtual const NormalVectorv_GetTraceNormal (const int edge) const final
 
virtual void v_ReOrientTracePhysMap (const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1)
 
virtual void v_TraceNormLen (const int traceid, NekDouble &h, NekDouble &p)
 
- 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)
 
virtual void v_MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_DivideByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_ComputeLaplacianMetric ()
 
Array< OneD, NekDoublev_GetMF (const int dir, const int shapedim, const StdRegions::VarCoeffMap &varcoeffs)
 
Array< OneD, NekDoublev_GetMFDiv (const int dir, const StdRegions::VarCoeffMap &varcoeffs)
 
Array< OneD, NekDoublev_GetMFMag (const int dir, const StdRegions::VarCoeffMap &varcoeffs)
 
virtual DNekMatSharedPtr v_BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
 
virtual DNekMatSharedPtr v_BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd)
 
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 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)
 
virtual void v_GetTraceQFactors (const int trace, Array< OneD, NekDouble > &outarray)
 
virtual void v_SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
virtual void v_SetUpPhysNormals (const int id)
 
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_GenTraceExp (const int traceid, ExpansionSharedPtr &exp)
 

Private Member Functions

 SegExp ()
 
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 51 of file SegExp.h.

Constructor & Destructor Documentation

◆ SegExp() [1/3]

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 59 of file SegExp.cpp.

61  : StdExpansion(Ba.GetNumModes(), 1, Ba),
62  StdExpansion1D(Ba.GetNumModes(), Ba), StdRegions::StdSegExp(Ba),
63  Expansion(geom), Expansion1D(geom),
65  std::bind(&SegExp::CreateMatrix, this, std::placeholders::_1),
66  std::string("SegExpMatrix")),
68  this, std::placeholders::_1),
69  std::string("SegExpStaticCondMatrix"))
70 {
71 }
Expansion1D(SpatialDomains::Geometry1DSharedPtr pGeom)
Definition: Expansion1D.h:60
DNekScalBlkMatSharedPtr CreateStaticCondMatrix(const MatrixKey &mkey)
Definition: Expansion.cpp:272
Expansion(SpatialDomains::GeometrySharedPtr pGeom)
Definition: Expansion.cpp:47
DNekScalMatSharedPtr CreateMatrix(const MatrixKey &mkey)
Definition: SegExp.cpp:1136
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: SegExp.h:240
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: SegExp.h:238
StdExpansion()
Default Constructor.

◆ SegExp() [2/3]

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

Copy Constructor.

Parameters
SExisting segment to duplicate.

Definition at line 77 of file SegExp.cpp.

78  : StdExpansion(S), StdExpansion1D(S), StdRegions::StdSegExp(S),
79  Expansion(S), Expansion1D(S), m_matrixManager(S.m_matrixManager),
80  m_staticCondMatrixManager(S.m_staticCondMatrixManager)
81 {
82 }

◆ ~SegExp()

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

Definition at line 87 of file SegExp.cpp.

88 {
89 }

◆ SegExp() [3/3]

Nektar::LocalRegions::SegExp::SegExp ( )
private

Member Function Documentation

◆ CreateMatrix()

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

Definition at line 1136 of file SegExp.cpp.

1137 {
1138  DNekScalMatSharedPtr returnval;
1139  NekDouble fac;
1141 
1143  "Geometric information is not set up");
1144 
1145  switch (mkey.GetMatrixType())
1146  {
1147  case StdRegions::eMass:
1148  {
1149  if ((m_metricinfo->GetGtype() == SpatialDomains::eDeformed) ||
1150  (mkey.GetNVarCoeff()))
1151  {
1152  fac = 1.0;
1153  goto UseLocRegionsMatrix;
1154  }
1155  else
1156  {
1157  fac = (m_metricinfo->GetJac(ptsKeys))[0];
1158  goto UseStdRegionsMatrix;
1159  }
1160  }
1161  break;
1162  case StdRegions::eInvMass:
1163  {
1164  if ((m_metricinfo->GetGtype() == SpatialDomains::eDeformed) ||
1165  (mkey.GetNVarCoeff()))
1166  {
1167  NekDouble one = 1.0;
1168  StdRegions::StdMatrixKey masskey(StdRegions::eMass,
1169  DetShapeType(), *this);
1170  DNekMatSharedPtr mat = GenMatrix(masskey);
1171  mat->Invert();
1172 
1173  returnval =
1175  }
1176  else
1177  {
1178  fac = 1.0 / (m_metricinfo->GetJac(ptsKeys))[0];
1179  goto UseStdRegionsMatrix;
1180  }
1181  }
1182  break;
1186  {
1187  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed ||
1188  mkey.GetNVarCoeff())
1189  {
1190  fac = 1.0;
1191  goto UseLocRegionsMatrix;
1192  }
1193  else
1194  {
1195  int dir = 0;
1196  switch (mkey.GetMatrixType())
1197  {
1199  dir = 0;
1200  break;
1202  ASSERTL1(m_geom->GetCoordim() >= 2,
1203  "Cannot call eWeakDeriv2 in a "
1204  "coordinate system which is not at "
1205  "least two-dimensional");
1206  dir = 1;
1207  break;
1209  ASSERTL1(m_geom->GetCoordim() == 3,
1210  "Cannot call eWeakDeriv2 in a "
1211  "coordinate system which is not "
1212  "three-dimensional");
1213  dir = 2;
1214  break;
1215  default:
1216  break;
1217  }
1218 
1219  MatrixKey deriv0key(StdRegions::eWeakDeriv0,
1220  mkey.GetShapeType(), *this);
1221 
1222  DNekMatSharedPtr WeakDerivStd = GetStdMatrix(deriv0key);
1223  fac = m_metricinfo->GetDerivFactors(ptsKeys)[dir][0] *
1224  m_metricinfo->GetJac(ptsKeys)[0];
1225 
1227  fac, WeakDerivStd);
1228  }
1229  }
1230  break;
1232  {
1233  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1234  {
1235  fac = 1.0;
1236  goto UseLocRegionsMatrix;
1237  }
1238  else
1239  {
1240  int coordim = m_geom->GetCoordim();
1241  fac = 0.0;
1242  for (int i = 0; i < coordim; ++i)
1243  {
1244  fac += m_metricinfo->GetDerivFactors(ptsKeys)[i][0] *
1245  m_metricinfo->GetDerivFactors(ptsKeys)[i][0];
1246  }
1247  fac *= m_metricinfo->GetJac(ptsKeys)[0];
1248  goto UseStdRegionsMatrix;
1249  }
1250  }
1251  break;
1253  {
1254  NekDouble factor = mkey.GetConstFactor(StdRegions::eFactorLambda);
1255  MatrixKey masskey(StdRegions::eMass, mkey.GetShapeType(), *this);
1256  DNekScalMat &MassMat = *(this->m_matrixManager[masskey]);
1257  MatrixKey lapkey(StdRegions::eLaplacian, mkey.GetShapeType(), *this,
1258  mkey.GetConstFactors(), mkey.GetVarCoeffs());
1259  DNekScalMat &LapMat = *(this->m_matrixManager[lapkey]);
1260 
1261  int rows = LapMat.GetRows();
1262  int cols = LapMat.GetColumns();
1263 
1264  DNekMatSharedPtr helm =
1266 
1267  NekDouble one = 1.0;
1268  (*helm) = LapMat + factor * MassMat;
1269 
1270  returnval =
1272  }
1273  break;
1278  {
1279  NekDouble one = 1.0;
1280 
1281  DNekMatSharedPtr mat = GenMatrix(mkey);
1282  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one, mat);
1283  }
1284  break;
1286  {
1287  NekDouble one = 1.0;
1288 
1289  // StdRegions::StdMatrixKey
1290  // hkey(StdRegions::eHybridDGHelmholtz,
1291  // DetShapeType(),*this,
1292  // mkey.GetConstant(0),
1293  // mkey.GetConstant(1));
1294  MatrixKey hkey(StdRegions::eHybridDGHelmholtz, DetShapeType(),
1295  *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
1296  DNekMatSharedPtr mat = GenMatrix(hkey);
1297 
1298  mat->Invert();
1299  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one, mat);
1300  }
1301  break;
1303  {
1304  DNekMatSharedPtr m_Ix;
1305  Array<OneD, NekDouble> coords(1, 0.0);
1306  StdRegions::ConstFactorMap factors = mkey.GetConstFactors();
1307  int vertex = (int)factors[StdRegions::eFactorGaussVertex];
1308 
1309  coords[0] = (vertex == 0) ? -1.0 : 1.0;
1310 
1311  m_Ix = m_base[0]->GetI(coords);
1312  returnval =
1314  }
1315  break;
1316 
1317  UseLocRegionsMatrix:
1318  {
1319  DNekMatSharedPtr mat = GenMatrix(mkey);
1320  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(fac, mat);
1321  }
1322  break;
1323  UseStdRegionsMatrix:
1324  {
1325  DNekMatSharedPtr mat = GetStdMatrix(mkey);
1326  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(fac, mat);
1327  }
1328  break;
1329  default:
1330  NEKERROR(ErrorUtil::efatal, "Matrix creation not defined");
1331  break;
1332  }
1333 
1334  return returnval;
1335 }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:272
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:272
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:273
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:611
const LibUtilities::PointsKeyVector GetPointsKeys() const
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:375
DNekMatSharedPtr GenMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:845
Array< OneD, LibUtilities::BasisSharedPtr > m_base
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:250
@ eNoGeomType
No type defined.
@ eDeformed
Geometry is curved or has non-constant factors.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:282
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::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 1408 of file SegExp.cpp.

1410 {
1411  // get Mass matrix inverse
1412  MatrixKey masskey(StdRegions::eInvMass, DetShapeType(), *this);
1413  DNekScalMatSharedPtr matsys = m_matrixManager[masskey];
1414 
1415  NekVector<NekDouble> in(m_ncoeffs, inarray, eCopy);
1416  NekVector<NekDouble> out(m_ncoeffs, outarray, eWrapper);
1417 
1418  out = (*matsys) * in;
1419 }

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 1367 of file SegExp.cpp.

1369 {
1370 
1371  int m;
1372  NekDouble sgn = 1;
1373 
1374  ASSERTL1(&inarray[0] != &outarray[0],
1375  "inarray and outarray can not be the same");
1376  switch (GetBasisType(0))
1377  {
1379  // Swap vertices
1380  outarray[0] = inarray[1];
1381  outarray[1] = inarray[0];
1382  // negate odd modes
1383  for (m = 2; m < m_ncoeffs; ++m)
1384  {
1385  outarray[m] = sgn * inarray[m];
1386  sgn = -sgn;
1387  }
1388  break;
1391  for (m = 0; m < m_ncoeffs; ++m)
1392  {
1393  outarray[m_ncoeffs - 1 - m] = inarray[m];
1394  }
1395  break;
1396  default:
1397  ASSERTL0(false, "This basis is not allowed in this method");
1398  break;
1399  }
1400 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:163
@ eGauss_Lagrange
Lagrange Polynomials using the Gauss points.
Definition: BasisType.h:59
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:58
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:50

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

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

Definition at line 673 of file SegExp.cpp.

675 {
676  size_t nquad = m_base[0]->GetNumPoints();
677 
679  {
680  switch (vertex)
681  {
682  case 0:
683  outarray[0] += inarray;
684  break;
685  case 1:
686  outarray[nquad - 1] += inarray;
687  break;
688  }
689  }
690  else
691  {
693  factors[StdRegions::eFactorGaussVertex] = vertex;
694 
695  StdRegions::StdMatrixKey key(StdRegions::eInterpGauss, DetShapeType(),
696  *this, factors);
697 
698  DNekScalMatSharedPtr mat_gauss = m_matrixManager[key];
699 
700  Vmath::Svtvp(nquad, inarray,
701  mat_gauss->GetOwnedMatrix()->GetPtr().get(), 1,
702  &outarray[0], 1, &outarray[0], 1);
703  }
704 }
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:213
@ eGaussGaussLegendre
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:48
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
Definition: Vmath.cpp:622

References Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::StdRegions::eFactorGaussVertex, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::StdRegions::eInterpGauss, Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::StdRegions::StdExpansion::m_base, m_matrixManager, and Vmath::Svtvp().

◆ v_ComputeTraceNormal()

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

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 854 of file SegExp.cpp.

855 {
856  int i;
857  const SpatialDomains::GeomFactorsSharedPtr &geomFactors =
858  GetGeom()->GetMetricInfo();
859  SpatialDomains::GeomType type = geomFactors->GetGtype();
860  const Array<TwoD, const NekDouble> &gmat =
861  geomFactors->GetDerivFactors(GetPointsKeys());
862  int nqe = 1;
863  int vCoordDim = GetCoordim();
864 
865  m_traceNormals[vertex] = Array<OneD, Array<OneD, NekDouble>>(vCoordDim);
866  Array<OneD, Array<OneD, NekDouble>> &normal = m_traceNormals[vertex];
867  for (i = 0; i < vCoordDim; ++i)
868  {
869  normal[i] = Array<OneD, NekDouble>(nqe);
870  }
871 
872  size_t nqb = nqe;
873  size_t nbnd = vertex;
874  m_elmtBndNormDirElmtLen[nbnd] = Array<OneD, NekDouble>{nqb, 0.0};
875  Array<OneD, NekDouble> &length = m_elmtBndNormDirElmtLen[nbnd];
876 
877  // Regular geometry case
878  if ((type == SpatialDomains::eRegular) ||
880  {
881  NekDouble vert;
882  // Set up normals
883  switch (vertex)
884  {
885  case 0:
886  for (i = 0; i < vCoordDim; ++i)
887  {
888  Vmath::Fill(nqe, -gmat[i][0], normal[i], 1);
889  }
890  break;
891  case 1:
892  for (i = 0; i < vCoordDim; ++i)
893  {
894  Vmath::Fill(nqe, gmat[i][0], normal[i], 1);
895  }
896  break;
897  default:
898  ASSERTL0(false, "point is out of range (point < 2)");
899  }
900 
901  // normalise
902  vert = 0.0;
903  for (i = 0; i < vCoordDim; ++i)
904  {
905  vert += normal[i][0] * normal[i][0];
906  }
907  vert = 1.0 / sqrt(vert);
908 
909  Vmath::Fill(nqb, vert, length, 1);
910 
911  for (i = 0; i < vCoordDim; ++i)
912  {
913  Vmath::Smul(nqe, vert, normal[i], 1, normal[i], 1);
914  }
915  }
916 }
std::map< int, NormalVector > m_traceNormals
Definition: Expansion.h:275
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:285
SpatialDomains::GeometrySharedPtr GetGeom() const
Definition: Expansion.cpp:166
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition: GeomFactors.h:62
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.cpp:248
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:291

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 1127 of file SegExp.cpp.

1128 {
1129  LibUtilities::BasisKey bkey = m_base[0]->GetBasisKey();
1132 
1133  return tmp->GetStdMatrix(mkey);
1134 }
std::shared_ptr< StdSegExp > StdSegExpSharedPtr
Definition: StdSegExp.h:46

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

◆ v_DropLocStaticCondMatrix()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1117 of file SegExp.cpp.

1118 {
1119  m_staticCondMatrixManager.DeleteObject(mkey);
1120 }

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 810 of file SegExp.cpp.

814 {
815  boost::ignore_unused(fromType);
816 
817  switch (m_base[0]->GetBasisType())
818  {
820  {
821  int fillorder = min((int)nummodes[mode_offset], m_ncoeffs);
822 
823  Vmath::Zero(m_ncoeffs, coeffs, 1);
824  Vmath::Vcopy(fillorder, &data[0], 1, &coeffs[0], 1);
825  }
826  break;
828  {
829  // Assume that input is also Gll_Lagrange
830  // but no way to check;
831  LibUtilities::PointsKey f0(nummodes[mode_offset],
833  LibUtilities::PointsKey t0(m_base[0]->GetNumModes(),
835  LibUtilities::Interp1D(f0, data, t0, coeffs);
836  }
837  break;
839  {
840  // Assume that input is also Gauss_Lagrange
841  // but no way to check;
842  LibUtilities::PointsKey f0(nummodes[mode_offset],
844  LibUtilities::PointsKey t0(m_base[0]->GetNumModes(),
846  LibUtilities::Interp1D(f0, data, t0, coeffs);
847  }
848  break;
849  default:
850  ASSERTL0(false, "basis is either not set up or not hierarchicial");
851  }
852 }
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:52
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:53
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255

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 357 of file SegExp.cpp.

359 {
360  if (m_base[0]->Collocation())
361  {
362  Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
363  }
364  else
365  {
366  v_IProductWRTBase(inarray, outarray);
367 
368  // get Mass matrix inverse
369  MatrixKey masskey(StdRegions::eInvMass, DetShapeType(), *this);
370  DNekScalMatSharedPtr matsys = m_matrixManager[masskey];
371 
372  // copy inarray in case inarray == outarray
373  NekVector<NekDouble> in(m_ncoeffs, outarray, eCopy);
374  NekVector<NekDouble> out(m_ncoeffs, outarray, eWrapper);
375 
376  out = (*matsys) * in;
377  }
378 }
virtual 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:483

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 380 of file SegExp.cpp.

383 {
384  if (m_base[0]->Collocation())
385  {
386  Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
387  }
388  else
389  {
390  int nInteriorDofs = m_ncoeffs - 2;
391  int offset = 0;
392 
393  switch (m_base[0]->GetBasisType())
394  {
396  {
397  offset = 1;
398  }
399  break;
401  {
402  nInteriorDofs = m_ncoeffs;
403  offset = 0;
404  }
405  break;
408  {
409  ASSERTL1(
410  m_base[0]->GetPointsType() ==
412  m_base[0]->GetPointsType() ==
414  "Cannot use FwdTrans_BndConstrained with these points.");
415  offset = 2;
416  }
417  break;
418  default:
419  ASSERTL0(false, "This type of FwdTrans is not defined"
420  "for this expansion type");
421  }
422 
423  fill(outarray.get(), outarray.get() + m_ncoeffs, 0.0);
424 
426  {
427 
428  outarray[GetVertexMap(0)] = inarray[0];
429  outarray[GetVertexMap(1)] = inarray[m_base[0]->GetNumPoints() - 1];
430 
431  if (m_ncoeffs > 2)
432  {
433  // ideally, we would like to have tmp0 to be replaced
434  // by outarray (currently MassMatrixOp does not allow
435  // aliasing)
436  Array<OneD, NekDouble> tmp0(m_ncoeffs);
437  Array<OneD, NekDouble> tmp1(m_ncoeffs);
438 
439  StdRegions::StdMatrixKey stdmasskey(StdRegions::eMass,
440  DetShapeType(), *this);
441  MassMatrixOp(outarray, tmp0, stdmasskey);
442  v_IProductWRTBase(inarray, tmp1);
443 
444  Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
445 
446  // get Mass matrix inverse (only of interior DOF)
447  MatrixKey masskey(StdRegions::eMass, DetShapeType(), *this);
448  DNekScalMatSharedPtr matsys =
449  (m_staticCondMatrixManager[masskey])->GetBlock(1, 1);
450 
451  Blas::Dgemv('N', nInteriorDofs, nInteriorDofs, matsys->Scale(),
452  &((matsys->GetOwnedMatrix())->GetPtr())[0],
453  nInteriorDofs, tmp1.get() + offset, 1, 0.0,
454  outarray.get() + offset, 1);
455  }
456  }
457  else
458  {
459  SegExp::v_FwdTrans(inarray, outarray);
460  }
461  }
462 }
virtual 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:357
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:760
int GetVertexMap(const int localVertexId, bool useCoeffPacking=false)
Definition: StdExpansion.h:687
static void Dgemv(const char &trans, const int &m, const int &n, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = A x where A[m x n].
Definition: Blas.hpp:246
@ ePolyEvenlySpaced
1D Evenly-spaced points using Lagrange polynomial
Definition: PointsType.h:75
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:51
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:419

References 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 1337 of file SegExp.cpp.

1338 {
1339  DNekMatSharedPtr returnval;
1340 
1341  switch (mkey.GetMatrixType())
1342  {
1349  returnval = Expansion1D::v_GenMatrix(mkey);
1350  break;
1351  default:
1352  returnval = StdSegExp::v_GenMatrix(mkey);
1353  break;
1354  }
1355 
1356  return returnval;
1357 }
virtual DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey)
Definition: Expansion1D.cpp:54

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

const LibUtilities::BasisSharedPtr & Nektar::LocalRegions::SegExp::v_GetBasis ( int  dir) const
protectedvirtual

Definition at line 793 of file SegExp.cpp.

794 {
795  return GetBasis(dir);
796 }
const LibUtilities::BasisSharedPtr & GetBasis(int dir) const
This function gets the shared point to basis in the dir direction.
Definition: StdExpansion.h:118

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

◆ 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 615 of file SegExp.cpp.

617 {
618  int i;
619 
620  ASSERTL1(Lcoords[0] >= -1.0 && Lcoords[0] <= 1.0,
621  "Local coordinates are not in region [-1,1]");
622 
623  m_geom->FillGeom();
624  for (i = 0; i < m_geom->GetCoordim(); ++i)
625  {
626  coords[i] = m_geom->GetCoord(i, Lcoords);
627  }
628 }

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

◆ v_GetCoordim()

int Nektar::LocalRegions::SegExp::v_GetCoordim ( void  )
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion1D.

Definition at line 767 of file SegExp.cpp.

768 {
769  return m_geom->GetCoordim();
770 }

References 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 630 of file SegExp.cpp.

633 {
634  Expansion::v_GetCoords(coords_0, coords_1, coords_2);
635 }
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3)
Definition: Expansion.cpp:524

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 759 of file SegExp.cpp.

760 {
761  LibUtilities::BasisKey bkey0(m_base[0]->GetBasisType(), 2,
762  m_base[0]->GetPointsKey());
763 
765 }

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 1122 of file SegExp.cpp.

1123 {
1124  return m_matrixManager[mkey];
1125 }

References m_matrixManager.

◆ v_GetLocStaticCondMatrix()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1112 of file SegExp.cpp.

1113 {
1114  return m_staticCondMatrixManager[mkey];
1115 }

References m_staticCondMatrixManager.

◆ v_GetNcoeffs()

int Nektar::LocalRegions::SegExp::v_GetNcoeffs ( void  ) const
protectedvirtual

Definition at line 788 of file SegExp.cpp.

789 {
790  return m_ncoeffs;
791 }

References Nektar::StdRegions::StdExpansion::m_ncoeffs.

◆ v_GetNumPoints()

int Nektar::LocalRegions::SegExp::v_GetNumPoints ( const int  dir) const
protectedvirtual

Definition at line 783 of file SegExp.cpp.

784 {
785  return GetNumPoints(dir);
786 }
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:226

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

◆ v_GetPhysNormals()

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

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 772 of file SegExp.cpp.

773 {
774  NEKERROR(ErrorUtil::efatal, "Got to SegExp");
775  return NullNekDouble1DArray;
776 }
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 719 of file SegExp.cpp.

720 {
721  int nquad = m_base[0]->GetNumPoints();
722 
723  ASSERTL1(vertex == 0 || vertex == 1, "Vertex value should be 0 or 1");
724 
725  map = Array<OneD, int>(1);
726 
727  map[0] = vertex == 0 ? 0 : nquad - 1;
728 }

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 706 of file SegExp.cpp.

710 {
711  boost::ignore_unused(EdgeExp, orient);
712 
713  NekDouble result;
714  v_GetVertexPhysVals(edge, inarray, result);
715  outarray[0] = result;
716 }
virtual void v_GetVertexPhysVals(const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray) override
Definition: SegExp.cpp:638

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 638 of file SegExp.cpp.

641 {
642  int nquad = m_base[0]->GetNumPoints();
643 
645  {
646  switch (vertex)
647  {
648  case 0:
649  outarray = inarray[0];
650  break;
651  case 1:
652  outarray = inarray[nquad - 1];
653  break;
654  }
655  }
656  else
657  {
659  factors[StdRegions::eFactorGaussVertex] = vertex;
660 
661  StdRegions::StdMatrixKey key(StdRegions::eInterpGauss, DetShapeType(),
662  *this, factors);
663 
664  DNekScalMatSharedPtr mat_gauss = m_matrixManager[key];
665 
666  outarray =
667  Blas::Ddot(nquad, mat_gauss->GetOwnedMatrix()->GetPtr().get(), 1,
668  &inarray[0], 1);
669  }
670 }
static double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = .
Definition: Blas.hpp:182

References Blas::Ddot(), Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::StdRegions::eFactorGaussVertex, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::StdRegions::eInterpGauss, 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 1013 of file SegExp.cpp.

1016 {
1017  int nquad = m_base[0]->GetNumPoints();
1018  const Array<TwoD, const NekDouble> &gmat =
1019  m_metricinfo->GetDerivFactors(GetPointsKeys());
1020  const NekDouble lambda = mkey.GetConstFactor(StdRegions::eFactorLambda);
1021 
1022  Array<OneD, NekDouble> physValues(nquad);
1023  Array<OneD, NekDouble> dPhysValuesdx(nquad);
1024  Array<OneD, NekDouble> wsp(m_ncoeffs);
1025 
1026  BwdTrans(inarray, physValues);
1027 
1028  // mass matrix operation
1029  v_IProductWRTBase((m_base[0]->GetBdata()), physValues, wsp, 1);
1030 
1031  // Laplacian matrix operation
1032  switch (m_geom->GetCoordim())
1033  {
1034  case 1:
1035  {
1036  PhysDeriv(physValues, dPhysValuesdx);
1037 
1038  // multiply with the proper geometric factors
1039  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1040  {
1041  Vmath::Vmul(nquad, &gmat[0][0], 1, dPhysValuesdx.get(), 1,
1042  dPhysValuesdx.get(), 1);
1043  }
1044  else
1045  {
1046  Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.get(), 1);
1047  }
1048  }
1049  break;
1050  case 2:
1051  {
1052  Array<OneD, NekDouble> dPhysValuesdy(nquad);
1053 
1054  PhysDeriv(physValues, dPhysValuesdx, dPhysValuesdy);
1055 
1056  // multiply with the proper geometric factors
1057  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1058  {
1059  Vmath::Vmul(nquad, &gmat[0][0], 1, dPhysValuesdx.get(), 1,
1060  dPhysValuesdx.get(), 1);
1061  Vmath::Vvtvp(nquad, &gmat[1][0], 1, dPhysValuesdy.get(), 1,
1062  dPhysValuesdx.get(), 1, dPhysValuesdx.get(), 1);
1063  }
1064  else
1065  {
1066  Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.get(), 1);
1067  Blas::Daxpy(nquad, gmat[1][0], dPhysValuesdy.get(), 1,
1068  dPhysValuesdx.get(), 1);
1069  }
1070  }
1071  break;
1072  case 3:
1073  {
1074  Array<OneD, NekDouble> dPhysValuesdy(nquad);
1075  Array<OneD, NekDouble> dPhysValuesdz(nquad);
1076 
1077  PhysDeriv(physValues, dPhysValuesdx, dPhysValuesdy, dPhysValuesdz);
1078 
1079  // multiply with the proper geometric factors
1080  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1081  {
1082  Vmath::Vmul(nquad, &gmat[0][0], 1, dPhysValuesdx.get(), 1,
1083  dPhysValuesdx.get(), 1);
1084  Vmath::Vvtvp(nquad, &gmat[1][0], 1, dPhysValuesdy.get(), 1,
1085  dPhysValuesdx.get(), 1, dPhysValuesdx.get(), 1);
1086  Vmath::Vvtvp(nquad, &gmat[2][0], 1, dPhysValuesdz.get(), 1,
1087  dPhysValuesdx.get(), 1, dPhysValuesdx.get(), 1);
1088  }
1089  else
1090  {
1091  Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.get(), 1);
1092  Blas::Daxpy(nquad, gmat[1][0], dPhysValuesdy.get(), 1,
1093  dPhysValuesdx.get(), 1);
1094  Blas::Daxpy(nquad, gmat[2][0], dPhysValuesdz.get(), 1,
1095  dPhysValuesdx.get(), 1);
1096  }
1097  }
1098  break;
1099  default:
1100  ASSERTL0(false, "Wrong number of dimensions");
1101  break;
1102  }
1103 
1104  v_IProductWRTBase(m_base[0]->GetDbdata(), dPhysValuesdx, outarray, 1);
1105  Blas::Daxpy(m_ncoeffs, lambda, wsp.get(), 1, outarray.get(), 1);
1106 }
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:432
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:850
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition: Blas.hpp:168
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
Definition: Blas.hpp:154
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:209
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:574

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 109 of file SegExp.cpp.

110 {
111  int nquad0 = m_base[0]->GetNumPoints();
112  Array<OneD, const NekDouble> jac = m_metricinfo->GetJac(GetPointsKeys());
113  NekDouble ival;
114  Array<OneD, NekDouble> tmp(nquad0);
115 
116  // multiply inarray with Jacobian
117  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
118  {
119  Vmath::Vmul(nquad0, jac, 1, inarray, 1, tmp, 1);
120  }
121  else
122  {
123  Vmath::Smul(nquad0, jac[0], inarray, 1, tmp, 1);
124  }
125 
126  // call StdSegExp version;
127  ival = StdSegExp::v_Integral(tmp);
128  // ival = StdSegExp::Integral(tmp);
129  return ival;
130 }

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 516 of file SegExp.cpp.

519 {
520  int nquad0 = m_base[0]->GetNumPoints();
521  Array<OneD, const NekDouble> jac = m_metricinfo->GetJac(GetPointsKeys());
522  Array<OneD, NekDouble> tmp(nquad0);
523 
524  // multiply inarray with Jacobian
525  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
526  {
527  Vmath::Vmul(nquad0, jac, 1, inarray, 1, tmp, 1);
528  }
529  else
530  {
531  Vmath::Smul(nquad0, jac[0], inarray, 1, tmp, 1);
532  }
533  StdSegExp::v_IProductWRTBase(base, tmp, outarray, coll_check);
534 }

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 483 of file SegExp.cpp.

485 {
486  v_IProductWRTBase(m_base[0]->GetBdata(), inarray, outarray, 1);
487 }

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

Referenced by v_FwdTrans(), v_FwdTransBndConstrained(), v_HelmholtzMatrixOp(), 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 536 of file SegExp.cpp.

539 {
540  ASSERTL1(dir < 3, "input dir is out of range");
541  ASSERTL1((dir == 2) ? m_geom->GetCoordim() == 3 : true,
542  "input dir is out of range");
543 
544  int nquad = m_base[0]->GetNumPoints();
545  const Array<TwoD, const NekDouble> &gmat =
546  m_metricinfo->GetDerivFactors(GetPointsKeys());
547 
548  Array<OneD, NekDouble> tmp1(nquad);
549 
550  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
551  {
552  Vmath::Vmul(nquad, gmat[dir], 1, inarray, 1, tmp1, 1);
553  }
554  else
555  {
556  Vmath::Smul(nquad, gmat[dir][0], inarray, 1, tmp1, 1);
557  }
558 
559  v_IProductWRTBase(m_base[0]->GetDbdata(), tmp1, outarray, 1);
560 }

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 922 of file SegExp.cpp.

925 {
926  boost::ignore_unused(mkey);
927 
928  int nquad = m_base[0]->GetNumPoints();
929  const Array<TwoD, const NekDouble> &gmat =
930  m_metricinfo->GetDerivFactors(GetPointsKeys());
931 
932  Array<OneD, NekDouble> physValues(nquad);
933  Array<OneD, NekDouble> dPhysValuesdx(nquad);
934 
935  BwdTrans(inarray, physValues);
936 
937  // Laplacian matrix operation
938  switch (m_geom->GetCoordim())
939  {
940  case 1:
941  {
942  PhysDeriv(physValues, dPhysValuesdx);
943 
944  // multiply with the proper geometric factors
945  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
946  {
947  Vmath::Vmul(nquad, &gmat[0][0], 1, dPhysValuesdx.get(), 1,
948  dPhysValuesdx.get(), 1);
949  }
950  else
951  {
952  Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.get(), 1);
953  }
954  }
955  break;
956  case 2:
957  {
958  Array<OneD, NekDouble> dPhysValuesdy(nquad);
959 
960  PhysDeriv(physValues, dPhysValuesdx, dPhysValuesdy);
961 
962  // multiply with the proper geometric factors
963  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
964  {
965  Vmath::Vmul(nquad, &gmat[0][0], 1, dPhysValuesdx.get(), 1,
966  dPhysValuesdx.get(), 1);
967  Vmath::Vvtvp(nquad, &gmat[1][0], 1, dPhysValuesdy.get(), 1,
968  dPhysValuesdx.get(), 1, dPhysValuesdx.get(), 1);
969  }
970  else
971  {
972  Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.get(), 1);
973  Blas::Daxpy(nquad, gmat[1][0], dPhysValuesdy.get(), 1,
974  dPhysValuesdx.get(), 1);
975  }
976  }
977  break;
978  case 3:
979  {
980  Array<OneD, NekDouble> dPhysValuesdy(nquad);
981  Array<OneD, NekDouble> dPhysValuesdz(nquad);
982 
983  PhysDeriv(physValues, dPhysValuesdx, dPhysValuesdy, dPhysValuesdz);
984 
985  // multiply with the proper geometric factors
986  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
987  {
988  Vmath::Vmul(nquad, &gmat[0][0], 1, dPhysValuesdx.get(), 1,
989  dPhysValuesdx.get(), 1);
990  Vmath::Vvtvp(nquad, &gmat[1][0], 1, dPhysValuesdy.get(), 1,
991  dPhysValuesdx.get(), 1, dPhysValuesdx.get(), 1);
992  Vmath::Vvtvp(nquad, &gmat[2][0], 1, dPhysValuesdz.get(), 1,
993  dPhysValuesdx.get(), 1, dPhysValuesdx.get(), 1);
994  }
995  else
996  {
997  Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.get(), 1);
998  Blas::Daxpy(nquad, gmat[1][0], dPhysValuesdy.get(), 1,
999  dPhysValuesdx.get(), 1);
1000  Blas::Daxpy(nquad, gmat[2][0], dPhysValuesdz.get(), 1,
1001  dPhysValuesdx.get(), 1);
1002  }
1003  }
1004  break;
1005  default:
1006  ASSERTL0(false, "Wrong number of dimensions");
1007  break;
1008  }
1009 
1010  v_IProductWRTBase(m_base[0]->GetDbdata(), dPhysValuesdx, outarray, 1);
1011 }

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

SpatialDomains::GeomType Nektar::LocalRegions::SegExp::v_MetricInfoType ( )
protectedvirtual

Definition at line 778 of file SegExp.cpp.

779 {
780  return m_metricinfo->GetGtype();
781 }

References Nektar::LocalRegions::Expansion::m_metricinfo.

◆ 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 580 of file SegExp.cpp.

583 {
584  NormVectorIProductWRTBase(Fvec[0], Fvec[1], outarray);
585 }
void NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:621

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 562 of file SegExp.cpp.

565 {
566  int nq = m_base[0]->GetNumPoints();
567  Array<OneD, NekDouble> Fn(nq);
568 
569  // @TODO: This routine no longer makes sense as a normal is not unique to an
570  // edge
571  const Array<OneD, const Array<OneD, NekDouble>> &normals =
572  GetLeftAdjacentElementExp()->GetTraceNormal(
574  Vmath::Vmul(nq, &Fx[0], 1, &normals[0][0], 1, &Fn[0], 1);
575  Vmath::Vvtvp(nq, &Fy[0], 1, &normals[1][0], 1, &Fn[0], 1, &Fn[0], 1);
576 
577  v_IProductWRTBase(Fn, outarray);
578 }
ExpansionSharedPtr GetLeftAdjacentElementExp() const
Definition: Expansion.h:433
int GetLeftAdjacentElementTrace() const
Definition: Expansion.h:446

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 798 of file SegExp.cpp.

799 {
800  return 2;
801 }

◆ v_NumDGBndryCoeffs()

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

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 803 of file SegExp.cpp.

804 {
805  return 2;
806 }

◆ 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 156 of file SegExp.cpp.

160 {
161  int nquad0 = m_base[0]->GetNumPoints();
162  Array<TwoD, const NekDouble> gmat =
163  m_metricinfo->GetDerivFactors(GetPointsKeys());
164  Array<OneD, NekDouble> diff(nquad0);
165 
166  // StdExpansion1D::PhysTensorDeriv(inarray,diff);
167  PhysTensorDeriv(inarray, diff);
168  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
169  {
170  if (out_d0.size())
171  {
172  Vmath::Vmul(nquad0, &gmat[0][0], 1, &diff[0], 1, &out_d0[0], 1);
173  }
174 
175  if (out_d1.size())
176  {
177  Vmath::Vmul(nquad0, &gmat[1][0], 1, &diff[0], 1, &out_d1[0], 1);
178  }
179 
180  if (out_d2.size())
181  {
182  Vmath::Vmul(nquad0, &gmat[2][0], 1, &diff[0], 1, &out_d2[0], 1);
183  }
184  }
185  else
186  {
187  if (out_d0.size())
188  {
189  Vmath::Smul(nquad0, gmat[0][0], diff, 1, out_d0, 1);
190  }
191 
192  if (out_d1.size())
193  {
194  Vmath::Smul(nquad0, gmat[1][0], diff, 1, out_d1, 1);
195  }
196 
197  if (out_d2.size())
198  {
199  Vmath::Smul(nquad0, gmat[2][0], diff, 1, out_d2, 1);
200  }
201  }
202 }
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 300 of file SegExp.cpp.

303 {
304  switch (dir)
305  {
306  case 0:
307  {
308  PhysDeriv(inarray, outarray, NullNekDouble1DArray,
310  }
311  break;
312  case 1:
313  {
314  PhysDeriv(inarray, NullNekDouble1DArray, outarray,
316  }
317  break;
318  case 2:
319  {
321  outarray);
322  }
323  break;
324  default:
325  {
326  ASSERTL1(false, "input dir is out of range");
327  }
328  break;
329  }
330 }

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 251 of file SegExp.cpp.

253 {
254  int nquad0 = m_base[0]->GetNumPoints();
255  Array<TwoD, const NekDouble> gmat =
256  m_metricinfo->GetDerivFactors(GetPointsKeys());
257  int coordim = m_geom->GetCoordim();
258  Array<OneD, NekDouble> out_dn_tmp(nquad0, 0.0);
259  switch (coordim)
260  {
261  case 2:
262 
263  Array<OneD, NekDouble> inarray_d0(nquad0);
264  Array<OneD, NekDouble> inarray_d1(nquad0);
265 
266  v_PhysDeriv(inarray, inarray_d0, inarray_d1);
267  Array<OneD, Array<OneD, NekDouble>> normals;
268  normals = Array<OneD, Array<OneD, NekDouble>>(coordim);
269  cout << "der_n" << endl;
270  for (int k = 0; k < coordim; ++k)
271  {
272  normals[k] = Array<OneD, NekDouble>(nquad0);
273  }
274  // @TODO: this routine no longer makes sense, since normals are not
275  // unique on
276  // an edge
277  // normals = GetMetricInfo()->GetNormal();
278  for (int i = 0; i < nquad0; i++)
279  {
280  cout << "nx= " << normals[0][i] << " ny=" << normals[1][i]
281  << endl;
282  }
284  "normal vectors do not exist: check if a"
285  "boundary region is defined as I ");
286  // \nabla u \cdot normal
287  Vmath::Vmul(nquad0, normals[0], 1, inarray_d0, 1, out_dn_tmp, 1);
288  Vmath::Vadd(nquad0, out_dn_tmp, 1, out_dn, 1, out_dn, 1);
289  Vmath::Zero(nquad0, out_dn_tmp, 1);
290  Vmath::Vmul(nquad0, normals[1], 1, inarray_d1, 1, out_dn_tmp, 1);
291  Vmath::Vadd(nquad0, out_dn_tmp, 1, out_dn, 1, out_dn, 1);
292 
293  for (int i = 0; i < nquad0; i++)
294  {
295  cout << "deps/dx =" << inarray_d0[i]
296  << " deps/dy=" << inarray_d1[i] << endl;
297  }
298  }
299 }
virtual 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:156
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.cpp:359

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 212 of file SegExp.cpp.

214 {
215  int nquad0 = m_base[0]->GetNumPoints();
216  int coordim = m_geom->GetCoordim();
217  Array<OneD, NekDouble> diff(nquad0);
218  // this operation is needed if you put out_ds==inarray
219  Vmath::Zero(nquad0, out_ds, 1);
220  switch (coordim)
221  {
222  case 2:
223  // diff= dU/de
224  Array<OneD, NekDouble> diff(nquad0);
225 
226  PhysTensorDeriv(inarray, diff);
227 
228  // get dS/de= (Jac)^-1
229  Array<OneD, NekDouble> Jac = m_metricinfo->GetJac(GetPointsKeys());
230  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
231  {
232  // calculate the derivative as (dU/de)*(Jac)^-1
233  Vmath::Vdiv(nquad0, diff, 1, Jac, 1, out_ds, 1);
234  }
235  else
236  {
237  NekDouble invJac = 1 / Jac[0];
238  Vmath::Smul(nquad0, invJac, diff, 1, out_ds, 1);
239  }
240  }
241 }
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.cpp:284

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

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 604 of file SegExp.cpp.

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

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 734 of file SegExp.cpp.

737 {
738 
739  if (dir == StdRegions::eBackwards)
740  {
741  if (&inarray[0] != &outarray[0])
742  {
743  Array<OneD, NekDouble> intmp(inarray);
744  ReverseCoeffsAndSign(intmp, outarray);
745  }
746  else
747  {
748  ReverseCoeffsAndSign(inarray, outarray);
749  }
750  }
751 }
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:1367

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 596 of file SegExp.cpp.

599 {
600  // Evaluate point in local (eta) coordinates.
601  return StdSegExp::v_PhysEvaluate(Lcoord, physvals);
602 }

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