Nektar++
Loading...
Searching...
No Matches
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, SpatialDomains::Geometry1D *geom)
 Constructor using BasisKey class for quadrature points and order definition.
 
 SegExp (const SegExp &S)
 Copy Constructor.
 
 ~SegExp () override=default
 
- Public Member Functions inherited from Nektar::StdRegions::StdSegExp
 StdSegExp (const LibUtilities::BasisKey &Ba)
 Constructor using BasisKey class for quadrature points and order definition.
 
 StdSegExp (const StdSegExp &T)=default
 
 ~StdSegExp () override=default
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion1D
 StdExpansion1D (int numcoeffs, const LibUtilities::BasisKey &Ba)
 
 StdExpansion1D ()=default
 
 StdExpansion1D (const StdExpansion1D &T)=default
 
 ~StdExpansion1D () override=default
 
void PhysTensorDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Evaluate the derivative \( d/d{\xi_1} \) at the physical quadrature points given by inarray and return in outarray.
 
NekDouble BaryTensorDeriv (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
 
NekDouble BaryTensorDeriv (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs, std::array< NekDouble, 6 > &secondOrderDerivs)
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion
 StdExpansion ()
 Default Constructor.
 
 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.
 
 StdExpansion (const StdExpansion &T)
 Copy Constructor.
 
virtual ~StdExpansion ()
 Destructor.
 
int GetNumBases () const
 This function returns the number of 1D bases used in the expansion.
 
const Array< OneD, const LibUtilities::BasisSharedPtr > & GetBase () const
 This function gets the shared point to basis.
 
const LibUtilities::BasisSharedPtrGetBasis (int dir) const
 This function gets the shared point to basis in the dir direction.
 
int GetNcoeffs (void) const
 This function returns the total number of coefficients used in the expansion.
 
int GetTotPoints () const
 This function returns the total number of quadrature points used in the element.
 
LibUtilities::BasisType GetBasisType (const int dir) const
 This function returns the type of basis used in the dir direction.
 
int GetBasisNumModes (const int dir) const
 This function returns the number of expansion modes in the dir direction.
 
int EvalBasisNumModesMax (void) const
 This function returns the maximum number of expansion modes over all local directions.
 
LibUtilities::PointsType GetPointsType (const int dir) const
 This function returns the type of quadrature points used in the dir direction.
 
int GetNumPoints (const int dir) const
 This function returns the number of quadrature points in the dir direction.
 
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.
 
int GetNverts () const
 This function returns the number of vertices of the expansion domain.
 
int GetTraceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th trace.
 
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.
 
const LibUtilities::BasisKey GetTraceBasisKey (const int i, int k=-1, bool UseGLL=false) const
 This function returns the basis key belonging to the i-th trace.
 
LibUtilities::PointsKey GetTracePointsKey (const int i, int k=-1) const
 This function returns the basis key belonging to the i-th trace.
 
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.
 
int GetNtraces () const
 Returns the number of trace elements connected to this element.
 
LibUtilities::ShapeType DetShapeType () const
 This function returns the shape of the expansion domain.
 
std::shared_ptr< StdExpansionGetStdExp () const
 
std::shared_ptr< StdExpansionGetLinStdExp (void) const
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion () const
 
bool IsNodalNonTensorialExp ()
 
void NodalToModal (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space.
 
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space.
 
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.
 
void FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 This function fills the array outarray with the mode-th mode of the expansion.
 
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
 
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.
 
void SetElmtId (const int id)
 Set the element id of this expansion when used in a list by returning value of m_elmt_id.
 
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
 
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
 
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}\)
 
void GeneralMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
void ExponentialFilter (Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
 
void LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
DNekMatSharedPtr GenMatrix (const StdMatrixKey &mkey)
 
void PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void PhysDeriv_s (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds)
 
void PhysDeriv_n (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn)
 
void PhysDirectionalDeriv (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &outarray)
 
void StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
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.
 
NekDouble PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
 This function evaluates the first derivative of the expansion at a single (arbitrary) point of the domain.
 
NekDouble PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs, std::array< NekDouble, 6 > &secondOrderDerivs)
 
NekDouble PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain.
 
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.
 
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta.
 
void LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi)
 Convert local collapsed coordinates eta into local cartesian coordinate xi.
 
void PhysInterp (std::shared_ptr< StdExpansion > fromExp, const Array< OneD, const NekDouble > &fromData, Array< OneD, NekDouble > &toData)
 interpolate from one set of quadrature points available from FromExp to the set of quadrature points in the current expansion. If the points are the same this routine will just copy the data
 
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.
 
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.
 
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.
 
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.
 
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.
 
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.
 
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::Geometry1D *pGeom)
 
 ~Expansion1D () override=default
 
void AddNormTraceInt (const int dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void AddHDGHelmholtzTraceTerms (const NekDouble tau, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
SpatialDomains::Geometry1DGetGeom1D () const
 
- Public Member Functions inherited from Nektar::LocalRegions::Expansion
 Expansion (SpatialDomains::Geometry *pGeom)
 
 Expansion (const Expansion &pSrc)
 
 ~Expansion () override
 
void SetTraceExp (const int traceid, ExpansionSharedPtr &f)
 
ExpansionSharedPtr GetTraceExp (const int traceid)
 
DNekScalMatSharedPtr GetLocMatrix (const LocalRegions::MatrixKey &mkey)
 
void DropLocMatrix (const LocalRegions::MatrixKey &mkey)
 
DNekScalMatSharedPtr GetLocMatrix (const StdRegions::MatrixType mtype, const StdRegions::ConstFactorMap &factors=StdRegions::NullConstFactorMap, const StdRegions::VarCoeffMap &varcoeffs=StdRegions::NullVarCoeffMap)
 
SpatialDomains::GeometryGetGeom () 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.
 
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).
 
void GetTracePhysVals (const int trace, const StdRegions::StdExpansionSharedPtr &TraceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=StdRegions::eNoOrientation)
 
void GetTracePhysMap (const int edge, Array< OneD, int > &outarray)
 
void ReOrientTracePhysMap (const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1)
 
const NormalVectorGetTraceNormal (const int id)
 
void ComputeTraceNormal (const int id)
 
const Array< OneD, const NekDouble > & GetPhysNormals (void)
 
void SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
void SetUpPhysNormals (const int trace)
 
void AddRobinMassMatrix (const int traceid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
void TraceNormLen (const int traceid, NekDouble &h, NekDouble &p)
 
void AddRobinTraceContribution (const int traceid, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs)
 
const Array< OneD, const NekDouble > & GetElmtBndNormDirElmtLen (const int nbnd) const
 
void StdDerivBaseOnTraceMat (Array< OneD, DNekMatSharedPtr > &DerivMat)
 
void PhysDerivBaseOnTraceMat (const int traceid, Array< OneD, DNekMatSharedPtr > &DerivMat)
 
void PhysBaseOnTraceMat (const int traceid, DNekMatSharedPtr &BdataMat)
 

Protected Member Functions

NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray) override
 Integrate the physical point list inarray over region and return the value.
 
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.
 
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.
 
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 \).
 
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 \).
 
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.
 
void v_FwdTransBndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Inner product of inarray over region with respect to the expansion basis (this)->_Base[0] and return in outarray.
 
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.
 
void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray) override
 
void v_NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray) override
 
NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals) override
 
NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &coord, const Array< OneD, const NekDouble > &physvals) override
 
NekDouble v_PhysEvalFirstDeriv (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
 
NekDouble v_PhysEvalFirstSecondDeriv (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs, std::array< NekDouble, 6 > &secondOrderDerivs) override
 
void v_GetCoord (const Array< OneD, const NekDouble > &Lcoords, Array< OneD, NekDouble > &coords) override
 
void v_GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
 
void v_GetVertexPhysVals (const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray) override
 
void v_GetTracePhysVals (const int edge, const StdRegions::StdExpansionSharedPtr &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient) override
 
void v_GetTracePhysMap (const int vertex, Array< OneD, int > &map) override
 
StdRegions::StdExpansionSharedPtr v_GetStdExp (void) const override
 
StdRegions::StdExpansionSharedPtr v_GetLinStdExp (void) const override
 
void v_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
int v_NumBndryCoeffs () const override
 
int v_NumDGBndryCoeffs () const override
 
void v_ComputeTraceNormal (const int vertex) override
 
void v_ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int mode_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType) override
 Unpack data from input file assuming it comes from.
 
const Array< OneD, const NekDouble > & v_GetPhysNormals () override
 
void v_LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
void v_LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
void v_HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
DNekMatSharedPtr v_GenMatrix (const StdRegions::StdMatrixKey &mkey) override
 
DNekScalMatSharedPtr CreateMatrix (const MatrixKey &mkey)
 
DNekMatSharedPtr v_CreateStdMatrix (const StdRegions::StdMatrixKey &mkey) override
 
DNekScalMatSharedPtr v_GetLocMatrix (const MatrixKey &mkey) override
 
void v_DropLocMatrix (const MatrixKey &mkey) override
 
DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const MatrixKey &mkey) override
 
void v_DropLocStaticCondMatrix (const MatrixKey &mkey) override
 
- Protected Member Functions inherited from Nektar::StdRegions::StdSegExp
NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray) override
 Integrate the physical point list inarray over region and return the value.
 
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.
 
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.
 
void v_StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray) override
 
void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Backward transform from coefficient space given in inarray and evaluate at the physical quadrature points outarray.
 
void v_BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_FwdTransBndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Forward transform from physical quadrature space stored in inarray and evaluate the expansion coefficients and store in outarray.
 
void v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Inner product of inarray over region with respect to the expansion basis (this)->m_base[0] and return in outarray.
 
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.
 
void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
 
void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta) override
 
void v_LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi) override
 
NekDouble v_PhysEvaluateBasis (const Array< OneD, const NekDouble > &coords, int mode) final
 
NekDouble v_PhysEvalFirstDeriv (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
 
NekDouble v_PhysEvalFirstSecondDeriv (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs, std::array< NekDouble, 6 > &secondOrderDerivs) override
 
void v_LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
 
void v_LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
 
void v_HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
 
void v_SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey) override
 
void v_ExponentialFilter (Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff) override
 
void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_FillMode (const int mode, Array< OneD, NekDouble > &outarray) override
 
void v_GetCoords (Array< OneD, NekDouble > &coords_0, Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2) override
 
void v_GetBoundaryMap (Array< OneD, unsigned int > &outarray) override
 
void v_GetInteriorMap (Array< OneD, unsigned int > &outarray) override
 
int v_GetVertexMap (int localVertexId, bool useCoeffPacking=false) override
 
int v_GetNverts () const final
 
int v_GetNtraces () const final
 
int v_GetTraceNcoeffs (const int i) const final
 
int v_GetTraceIntNcoeffs (const int i) const final
 
int v_GetTraceNumPoints (const int i) const final
 
int v_NumBndryCoeffs () const override
 
int v_NumDGBndryCoeffs () const override
 
bool v_IsBoundaryInteriorExpansion () const override
 
int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset) override
 
LibUtilities::ShapeType v_DetShapeType () const override
 Return Shape of region, using ShapeType enum list. i.e. Segment.
 
DNekMatSharedPtr v_GenMatrix (const StdMatrixKey &mkey) override
 
DNekMatSharedPtr v_CreateStdMatrix (const StdMatrixKey &mkey) override
 
void v_GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true) override
 
void v_ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_GetTraceCoeffMap (const unsigned int traceid, Array< OneD, unsigned int > &maparray) override
 Get the map of the coefficient location to teh local trace coefficients.
 
void v_GetElmtTraceToTraceMap (const unsigned int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient, int P, int Q) override
 
void v_GetTraceToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient, int P, int Q) override
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion1D
NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals) override
 
void v_PhysInterp (std::shared_ptr< StdExpansion > fromExp, const Array< OneD, const NekDouble > &fromData, Array< OneD, NekDouble > &toData) override
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion
DNekMatSharedPtr CreateStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr CreateStdStaticCondMatrix (const StdMatrixKey &mkey)
 Create the static condensation of a matrix when using a boundary interior decomposition.
 
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase_SumFac (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GeneralMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
 
void LaplacianMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp_MatFree (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void HelmholtzMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals, NekDouble &deriv, NekDouble &deriv2)
 This function performs the barycentric interpolation of the polynomial stored in coord at a point physvals using barycentric interpolation weights in direction.
 
template<int DIR>
NekDouble BaryEvaluateBasis (const NekDouble &coord, const int &mode)
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals)
 Helper function to pass an unused value by reference into BaryEvaluate.
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals, NekDouble &deriv)
 
- Protected Member Functions inherited from Nektar::LocalRegions::Expansion1D
void v_AddRobinMassMatrix (const int vert, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat) override
 
void v_AddRobinTraceContribution (const int vert, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs) override
 
NekDouble v_VectorFlux (const Array< OneD, Array< OneD, NekDouble > > &vec) override
 
void v_NormalTraceDerivFactors (Array< OneD, Array< OneD, NekDouble > > &factors, Array< OneD, Array< OneD, NekDouble > > &d0factors, Array< OneD, Array< OneD, NekDouble > > &d1factors) override
 : This method gets all of the factors which are required as part of the Gradient Jump Penalty stabilisation and involves the product of the normal and geometric factors along the element trace.
 
void v_ReOrientTracePhysMap (const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1) override
 
void v_TraceNormLen (const int traceid, NekDouble &h, NekDouble &p) override
 
- Protected Member Functions inherited from Nektar::LocalRegions::Expansion
void ComputeLaplacianMetric ()
 
void ComputeQuadratureMetric ()
 
void ComputeGmatcdotMF (const Array< TwoD, const NekDouble > &df, const Array< OneD, const NekDouble > &direction, Array< OneD, Array< OneD, NekDouble > > &dfdir)
 
Array< OneD, NekDoubleGetMF (const int dir, const int shapedim, const StdRegions::VarCoeffMap &varcoeffs)
 
Array< OneD, NekDoubleGetMFDiv (const int dir, const StdRegions::VarCoeffMap &varcoeffs)
 
Array< OneD, NekDoubleGetMFMag (const int dir, const StdRegions::VarCoeffMap &varcoeffs)
 
void v_MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
virtual void v_DivideByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_ComputeLaplacianMetric ()
 
int v_GetCoordim () const override
 
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_GenTraceExp (const int traceid, ExpansionSharedPtr &exp)
 

Private Member Functions

void ReverseCoeffsAndSign (const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Reverse the coefficients in a boundary interior expansion this routine is of use when we need the segment coefficients corresponding to a expansion in the reverse coordinate direction.
 
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::Geometrym_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)
 

Detailed Description

Defines a Segment local expansion.

Definition at line 49 of file SegExp.h.

Constructor & Destructor Documentation

◆ SegExp() [1/2]

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

Constructor using BasisKey class for quadrature points and order definition.

Parameters
BaBasis key of segment expansion.
geomDescription of geometry.

Definition at line 55 of file SegExp.cpp.

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

◆ SegExp() [2/2]

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

Copy Constructor.

Parameters
SExisting segment to duplicate.

Definition at line 73 of file SegExp.cpp.

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

◆ ~SegExp()

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

Member Function Documentation

◆ CreateMatrix()

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

Definition at line 1127 of file SegExp.cpp.

1128{
1129 DNekScalMatSharedPtr returnval;
1130 NekDouble fac;
1132
1134 "Geometric information is not set up");
1135
1136 switch (mkey.GetMatrixType())
1137 {
1138 case StdRegions::eMass:
1139 {
1140 if ((m_metricinfo->GetGtype() == SpatialDomains::eDeformed) ||
1141 (mkey.GetNVarCoeff()))
1142 {
1143 fac = 1.0;
1144 goto UseLocRegionsMatrix;
1145 }
1146 else
1147 {
1148 fac = (m_metricinfo->GetJac(ptsKeys))[0];
1149 goto UseStdRegionsMatrix;
1150 }
1151 }
1152 break;
1154 {
1155 if ((m_metricinfo->GetGtype() == SpatialDomains::eDeformed) ||
1156 (mkey.GetNVarCoeff()))
1157 {
1158 NekDouble one = 1.0;
1159 StdRegions::StdMatrixKey masskey(StdRegions::eMass,
1160 DetShapeType(), *this);
1161 DNekMatSharedPtr mat = GenMatrix(masskey);
1162 mat->Invert();
1163
1164 returnval =
1166 }
1167 else
1168 {
1169 fac = 1.0 / (m_metricinfo->GetJac(ptsKeys))[0];
1170 goto UseStdRegionsMatrix;
1171 }
1172 }
1173 break;
1177 {
1178 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed ||
1179 mkey.GetNVarCoeff())
1180 {
1181 fac = 1.0;
1182 goto UseLocRegionsMatrix;
1183 }
1184 else
1185 {
1186 int dir = 0;
1187 switch (mkey.GetMatrixType())
1188 {
1190 dir = 0;
1191 break;
1193 ASSERTL1(m_geom->GetCoordim() >= 2,
1194 "Cannot call eWeakDeriv2 in a "
1195 "coordinate system which is not at "
1196 "least two-dimensional");
1197 dir = 1;
1198 break;
1200 ASSERTL1(m_geom->GetCoordim() == 3,
1201 "Cannot call eWeakDeriv2 in a "
1202 "coordinate system which is not "
1203 "three-dimensional");
1204 dir = 2;
1205 break;
1206 default:
1207 break;
1208 }
1209
1210 MatrixKey deriv0key(StdRegions::eWeakDeriv0,
1211 mkey.GetShapeType(), *this);
1212
1213 DNekMatSharedPtr WeakDerivStd = GetStdMatrix(deriv0key);
1214 fac = m_metricinfo->GetDerivFactors(ptsKeys)[dir][0] *
1215 m_metricinfo->GetJac(ptsKeys)[0];
1216
1218 fac, WeakDerivStd);
1219 }
1220 }
1221 break;
1223 {
1224 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1225 {
1226 fac = 1.0;
1227 goto UseLocRegionsMatrix;
1228 }
1229 else
1230 {
1231 int coordim = m_geom->GetCoordim();
1232 fac = 0.0;
1233 for (int i = 0; i < coordim; ++i)
1234 {
1235 fac += m_metricinfo->GetDerivFactors(ptsKeys)[i][0] *
1236 m_metricinfo->GetDerivFactors(ptsKeys)[i][0];
1237 }
1238 fac *= m_metricinfo->GetJac(ptsKeys)[0];
1239 goto UseStdRegionsMatrix;
1240 }
1241 }
1242 break;
1244 {
1245 if ((m_metricinfo->GetGtype() == SpatialDomains::eDeformed) ||
1246 (mkey.GetNVarCoeff()))
1247 {
1248 fac = 1.0;
1249 goto UseLocRegionsMatrix;
1250 }
1251 else
1252 {
1253 fac = (m_metricinfo->GetJac(ptsKeys))[0];
1254 goto UseStdRegionsMatrix;
1255 }
1256 }
1257 break;
1259 {
1260 NekDouble factor = mkey.GetConstFactor(StdRegions::eFactorLambda);
1261 MatrixKey masskey(StdRegions::eMass, mkey.GetShapeType(), *this);
1262 DNekScalMat &MassMat = *(this->m_matrixManager[masskey]);
1263 MatrixKey lapkey(StdRegions::eLaplacian, mkey.GetShapeType(), *this,
1264 mkey.GetConstFactors(), mkey.GetVarCoeffs());
1265 DNekScalMat &LapMat = *(this->m_matrixManager[lapkey]);
1266
1267 int rows = LapMat.GetRows();
1268 int cols = LapMat.GetColumns();
1269
1270 DNekMatSharedPtr helm =
1272
1273 NekDouble one = 1.0;
1274 (*helm) = LapMat + factor * MassMat;
1275
1276 returnval =
1278 }
1279 break;
1281 {
1282 NekDouble lambda = mkey.GetConstFactor(StdRegions::eFactorLambda);
1283
1284 // Construct mass matrix
1285 // Check for mass-specific varcoeffs to avoid unncessary
1286 // re-computation of the elemental matrix every time step
1288 if (mkey.HasVarCoeff(StdRegions::eVarCoeffMass))
1289 {
1290 massVarcoeffs[StdRegions::eVarCoeffMass] =
1291 mkey.GetVarCoeff(StdRegions::eVarCoeffMass);
1292 }
1293 MatrixKey masskey(StdRegions::eMass, mkey.GetShapeType(), *this,
1294 mkey.GetConstFactors(), massVarcoeffs);
1295 DNekScalMat &MassMat = *GetLocMatrix(masskey);
1296
1297 // Construct advection matrix
1298 // Check for varcoeffs not required;
1299 // assume advection velocity is always time-dependent
1300 MatrixKey advkey(mkey, StdRegions::eLinearAdvection);
1301 DNekScalMat &AdvMat = *GetLocMatrix(advkey);
1302
1303 int rows = MassMat.GetRows();
1304 int cols = MassMat.GetColumns();
1305
1306 DNekMatSharedPtr adr =
1308
1309 NekDouble one = 1.0;
1310 (*adr) = -lambda * MassMat + AdvMat;
1311
1313
1314 // Clear memory for time-dependent matrices
1315 DropLocMatrix(advkey);
1316 if (!massVarcoeffs.empty())
1317 {
1318 DropLocMatrix(masskey);
1319 }
1320 }
1321 break;
1323 {
1324 NekDouble lambda = mkey.GetConstFactor(StdRegions::eFactorLambda);
1325
1326 // Construct mass matrix
1327 // Check for mass-specific varcoeffs to avoid unncessary
1328 // re-computation of the elemental matrix every time step
1330 if (mkey.HasVarCoeff(StdRegions::eVarCoeffMass))
1331 {
1332 massVarcoeffs[StdRegions::eVarCoeffMass] =
1333 mkey.GetVarCoeff(StdRegions::eVarCoeffMass);
1334 }
1335 MatrixKey masskey(StdRegions::eMass, mkey.GetShapeType(), *this,
1336 mkey.GetConstFactors(), massVarcoeffs);
1337 DNekScalMat &MassMat = *GetLocMatrix(masskey);
1338
1339 // Construct laplacian matrix (Check for varcoeffs)
1340 // Take all varcoeffs if one or more are detected
1341 // TODO We might want to have a map
1342 // from MatrixType to Vector of Varcoeffs and vice-versa
1344 if ((mkey.HasVarCoeff(StdRegions::eVarCoeffLaplacian)) ||
1345 (mkey.HasVarCoeff(StdRegions::eVarCoeffD00)) ||
1346 (mkey.HasVarCoeff(StdRegions::eVarCoeffD01)) ||
1347 (mkey.HasVarCoeff(StdRegions::eVarCoeffD10)) ||
1348 (mkey.HasVarCoeff(StdRegions::eVarCoeffD02)) ||
1349 (mkey.HasVarCoeff(StdRegions::eVarCoeffD20)) ||
1350 (mkey.HasVarCoeff(StdRegions::eVarCoeffD11)) ||
1351 (mkey.HasVarCoeff(StdRegions::eVarCoeffD12)) ||
1352 (mkey.HasVarCoeff(StdRegions::eVarCoeffD21)) ||
1353 (mkey.HasVarCoeff(StdRegions::eVarCoeffD22)))
1354 {
1355 lapVarcoeffs = mkey.GetVarCoeffs();
1356 }
1357 MatrixKey lapkey(StdRegions::eLaplacian, mkey.GetShapeType(), *this,
1358 mkey.GetConstFactors(), lapVarcoeffs);
1359 DNekScalMat &LapMat = *GetLocMatrix(lapkey);
1360
1361 // Construct advection matrix
1362 // Check for varcoeffs not required;
1363 // assume advection velocity is always time-dependent
1364 MatrixKey advkey(mkey, StdRegions::eLinearAdvection);
1365 DNekScalMat &AdvMat = *GetLocMatrix(advkey);
1366
1367 int rows = LapMat.GetRows();
1368 int cols = LapMat.GetColumns();
1369
1370 DNekMatSharedPtr adr =
1372
1373 NekDouble one = 1.0;
1374 (*adr) = LapMat - lambda * MassMat + AdvMat;
1375
1377
1378 // Clear memory for time-dependent matrices
1379 DropLocMatrix(advkey);
1380 if (!massVarcoeffs.empty())
1381 {
1382 DropLocMatrix(masskey);
1383 }
1384 if (!lapVarcoeffs.empty())
1385 {
1386 DropLocMatrix(lapkey);
1387 }
1388 }
1389 break;
1394 {
1395 NekDouble one = 1.0;
1396
1397 DNekMatSharedPtr mat = GenMatrix(mkey);
1399 }
1400 break;
1402 {
1403 NekDouble one = 1.0;
1404
1405 // StdRegions::StdMatrixKey
1406 // hkey(StdRegions::eHybridDGHelmholtz,
1407 // DetShapeType(),*this,
1408 // mkey.GetConstant(0),
1409 // mkey.GetConstant(1));
1411 *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
1412 DNekMatSharedPtr mat = GenMatrix(hkey);
1413
1414 mat->Invert();
1416 }
1417 break;
1419 {
1420 DNekMatSharedPtr m_Ix;
1421 Array<OneD, NekDouble> coords(1, 0.0);
1422 StdRegions::ConstFactorMap factors = mkey.GetConstFactors();
1423 int vertex = (int)factors[StdRegions::eFactorGaussVertex];
1424
1425 coords[0] = (vertex == 0) ? -1.0 : 1.0;
1426
1427 m_Ix = m_base[0]->GetI(coords);
1428 returnval =
1430 }
1431 break;
1432
1433 UseLocRegionsMatrix:
1434 {
1435 DNekMatSharedPtr mat = GenMatrix(mkey);
1437 }
1438 break;
1439 UseStdRegionsMatrix:
1440 {
1441 DNekMatSharedPtr mat = GetStdMatrix(mkey);
1443 }
1444 break;
1445 default:
1446 {
1447 NekDouble one = 1.0;
1448 DNekMatSharedPtr mat = GenMatrix(mkey);
1449
1451 }
1452 break;
1453 }
1454
1455 return returnval;
1456}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
void DropLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition Expansion.cpp:90
SpatialDomains::Geometry * m_geom
Definition Expansion.h:279
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition Expansion.h:280
DNekScalMatSharedPtr GetLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition Expansion.cpp:84
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
int GetCoordim() const
Return the coordinate dimension of this object (i.e. the dimension of the space in which this object ...
Definition Geometry.h:279
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
const LibUtilities::PointsKeyVector GetPointsKeys() const
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
DNekMatSharedPtr GenMatrix(const StdMatrixKey &mkey)
Array< OneD, LibUtilities::BasisSharedPtr > m_base
std::vector< PointsKey > PointsKeyVector
Definition Points.h:231
@ eNoGeomType
No type defined.
@ eDeformed
Geometry is curved or has non-constant factors.
std::map< ConstFactorType, NekDouble > ConstFactorMap
static VarCoeffMap NullVarCoeffMap
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
StdRegions::ConstFactorMap factors
NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag > DNekScalMat
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekMat > DNekMatSharedPtr

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL1, ASSERTL2, Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::LocalRegions::Expansion::DropLocMatrix(), Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::eFactorGaussVertex, Nektar::StdRegions::eFactorLambda, 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::eLinearAdvection, Nektar::StdRegions::eLinearAdvectionDiffusionReaction, Nektar::StdRegions::eLinearAdvectionReaction, Nektar::StdRegions::eMass, Nektar::SpatialDomains::eNoGeomType, Nektar::StdRegions::eVarCoeffD00, Nektar::StdRegions::eVarCoeffD01, Nektar::StdRegions::eVarCoeffD02, Nektar::StdRegions::eVarCoeffD10, Nektar::StdRegions::eVarCoeffD11, Nektar::StdRegions::eVarCoeffD12, Nektar::StdRegions::eVarCoeffD20, Nektar::StdRegions::eVarCoeffD21, Nektar::StdRegions::eVarCoeffD22, Nektar::StdRegions::eVarCoeffLaplacian, Nektar::StdRegions::eVarCoeffMass, Nektar::StdRegions::eWeakDeriv0, Nektar::StdRegions::eWeakDeriv1, Nektar::StdRegions::eWeakDeriv2, Nektar::StdRegions::StdExpansion::GenMatrix(), Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdMatrixKey::GetConstFactors(), Nektar::SpatialDomains::Geometry::GetCoordim(), Nektar::LocalRegions::Expansion::GetLocMatrix(), Nektar::StdRegions::StdMatrixKey::GetMatrixType(), Nektar::StdRegions::StdMatrixKey::GetNVarCoeff(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdMatrixKey::GetShapeType(), Nektar::StdRegions::StdExpansion::GetStdMatrix(), Nektar::StdRegions::StdMatrixKey::GetVarCoeff(), Nektar::StdRegions::StdMatrixKey::GetVarCoeffs(), Nektar::StdRegions::StdMatrixKey::HasVarCoeff(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_geom, m_matrixManager, Nektar::LocalRegions::Expansion::m_metricinfo, and Nektar::StdRegions::NullVarCoeffMap.

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

1531{
1532 // get Mass matrix inverse
1533 MatrixKey masskey(StdRegions::eInvMass, DetShapeType(), *this);
1534 DNekScalMatSharedPtr matsys = m_matrixManager[masskey];
1535
1536 NekVector<NekDouble> in(m_ncoeffs, inarray, eCopy);
1537 NekVector<NekDouble> out(m_ncoeffs, outarray, eWrapper);
1538
1539 out = (*matsys) * in;
1540}

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

1490{
1491
1492 int m;
1493 NekDouble sgn = 1;
1494
1495 ASSERTL1(&inarray[0] != &outarray[0],
1496 "inarray and outarray can not be the same");
1497 switch (GetBasisType(0))
1498 {
1500 // Swap vertices
1501 outarray[0] = inarray[1];
1502 outarray[1] = inarray[0];
1503 // negate odd modes
1504 for (m = 2; m < m_ncoeffs; ++m)
1505 {
1506 outarray[m] = sgn * inarray[m];
1507 sgn = -sgn;
1508 }
1509 break;
1512 for (m = 0; m < m_ncoeffs; ++m)
1513 {
1514 outarray[m_ncoeffs - 1 - m] = inarray[m];
1515 }
1516 break;
1517 default:
1518 ASSERTL0(false, "This basis is not allowed in this method");
1519 break;
1520 }
1521}
#define ASSERTL0(condition, msg)
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
@ eGauss_Lagrange
Lagrange Polynomials using the Gauss points.
Definition BasisType.h:57
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition BasisType.h:56
@ eModified_A
Principle Modified Functions .
Definition BasisType.h:48

References ASSERTL0, ASSERTL1, Nektar::LibUtilities::eGauss_Lagrange, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), and Nektar::StdRegions::StdExpansion::m_ncoeffs.

Referenced by v_SetCoeffsToOrientation().

◆ v_ComputeTraceNormal()

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

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 833 of file SegExp.cpp.

834{
835 int i;
836 const SpatialDomains::GeomFactorsSharedPtr &geomFactors =
838 SpatialDomains::GeomType type = geomFactors->GetGtype();
839 const Array<TwoD, const NekDouble> &gmat =
840 geomFactors->GetDerivFactors(GetPointsKeys());
841 int nqe = 1;
842 int vCoordDim = GetCoordim();
843
844 m_traceNormals[vertex] = Array<OneD, Array<OneD, NekDouble>>(vCoordDim);
845 Array<OneD, Array<OneD, NekDouble>> &normal = m_traceNormals[vertex];
846 for (i = 0; i < vCoordDim; ++i)
847 {
848 normal[i] = Array<OneD, NekDouble>(nqe);
849 }
850
851 size_t nqb = nqe;
852 size_t nbnd = vertex;
853 m_elmtBndNormDirElmtLen[nbnd] = Array<OneD, NekDouble>{nqb, 0.0};
854 Array<OneD, NekDouble> &length = m_elmtBndNormDirElmtLen[nbnd];
855
856 // Regular geometry case
857 if ((type == SpatialDomains::eRegular) ||
859 {
860 NekDouble vert;
861 // Set up normals
862 switch (vertex)
863 {
864 case 0:
865 for (i = 0; i < vCoordDim; ++i)
866 {
867 Vmath::Fill(nqe, -gmat[i][0], normal[i], 1);
868 }
869 break;
870 case 1:
871 for (i = 0; i < vCoordDim; ++i)
872 {
873 Vmath::Fill(nqe, gmat[i][0], normal[i], 1);
874 }
875 break;
876 default:
877 ASSERTL0(false, "point is out of range (point < 2)");
878 }
879
880 // normalise
881 vert = 0.0;
882 for (i = 0; i < vCoordDim; ++i)
883 {
884 vert += normal[i][0] * normal[i][0];
885 }
886 vert = 1.0 / sqrt(vert);
887
888 Vmath::Fill(nqb, vert, length, 1);
889
890 for (i = 0; i < vCoordDim; ++i)
891 {
892 Vmath::Smul(nqe, vert, normal[i], 1, normal[i], 1);
893 }
894 }
895}
std::map< int, NormalVector > m_traceNormals
Definition Expansion.h:282
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:292
SpatialDomains::Geometry * GetGeom() const
GeomFactorsSharedPtr GetMetricInfo()
Get the geometric factors for this object.
Definition Geometry.h:306
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition GeomFactors.h:58
GeomType
Indicates the type of element geometry.
@ eRegular
Geometry is straight-sided with constant geometric factors.
@ eMovingRegular
Currently unused.
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition Vmath.hpp:100
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition Vmath.hpp:54
scalarT< T > sqrt(scalarT< T > in)
Definition scalar.hpp:290

References ASSERTL0, Nektar::SpatialDomains::eMovingRegular, Nektar::SpatialDomains::eRegular, Vmath::Fill(), Nektar::StdRegions::StdExpansion::GetCoordim(), Nektar::LocalRegions::Expansion::GetGeom(), Nektar::SpatialDomains::Geometry::GetMetricInfo(), 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::StdExpansion.

Definition at line 1118 of file SegExp.cpp.

1119{
1120 LibUtilities::BasisKey bkey = m_base[0]->GetBasisKey();
1123
1124 return tmp->GetStdMatrix(mkey);
1125}
std::shared_ptr< StdSegExp > StdSegExpSharedPtr
Definition StdSegExp.h:214

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

◆ v_DropLocMatrix()

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

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1113 of file SegExp.cpp.

1114{
1115 m_matrixManager.DeleteObject(mkey);
1116}

References m_matrixManager.

◆ v_DropLocStaticCondMatrix()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1103 of file SegExp.cpp.

1104{
1105 m_staticCondMatrixManager.DeleteObject(mkey);
1106}

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

795{
796 switch (m_base[0]->GetBasisType())
797 {
799 {
800 int fillorder = min((int)nummodes[mode_offset], m_ncoeffs);
801
802 Vmath::Zero(m_ncoeffs, coeffs, 1);
803 Vmath::Vcopy(fillorder, &data[0], 1, &coeffs[0], 1);
804 }
805 break;
807 {
808 // Assume that input is also Gll_Lagrange
809 // but no way to check;
810 LibUtilities::PointsKey f0(nummodes[mode_offset],
812 LibUtilities::PointsKey t0(m_base[0]->GetNumModes(),
814 LibUtilities::Interp1D(f0, data, t0, coeffs);
815 }
816 break;
818 {
819 // Assume that input is also Gauss_Lagrange
820 // but no way to check;
821 LibUtilities::PointsKey f0(nummodes[mode_offset],
823 LibUtilities::PointsKey t0(m_base[0]->GetNumModes(),
825 LibUtilities::Interp1D(f0, data, t0, coeffs);
826 }
827 break;
828 default:
829 ASSERTL0(false, "basis is either not set up or not hierarchicial");
830 }
831}
void Interp1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
this function interpolates a 1D function evaluated at the quadrature points of the basis fbasis0 to ...
Definition Interp.cpp:47
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition PointsType.h:51
@ eGaussGaussLegendre
1D Gauss-Gauss-Legendre quadrature points
Definition PointsType.h:46
void Zero(int n, T *x, const int incx)
Zero vector.
Definition Vmath.hpp:273
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition Vmath.hpp:825
scalarT< T > min(scalarT< T > lhs, scalarT< T > rhs)
Definition scalar.hpp:300

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, tinysimd::min(), 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.

Implements Nektar::StdRegions::StdExpansion.

Definition at line 346 of file SegExp.cpp.

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

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::StdExpansion.

Definition at line 369 of file SegExp.cpp.

372{
373 if (m_base[0]->Collocation())
374 {
375 Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
376 }
377 else
378 {
379 int nInteriorDofs = m_ncoeffs - 2;
380 int offset = 0;
381 bool hasEndPoints = true;
382 bool hasEndModes = true;
383
384 switch (m_base[0]->GetBasisType())
385 {
387 {
388 nInteriorDofs = m_ncoeffs - 2;
389 offset = 1;
390 hasEndModes = true;
391 }
392 break;
395 {
396 nInteriorDofs = m_ncoeffs;
397 offset = 0;
398 hasEndModes = false;
399 }
400 break;
403 {
404 nInteriorDofs = m_ncoeffs - 2;
405 offset = 2;
406 hasEndModes = true;
407 }
408 break;
409 default:
410 ASSERTL0(false, "This type of FwdTrans is not defined"
411 "for this expansion type");
412 }
413
414 switch (m_base[0]->GetPointsType())
415 {
418 case LibUtilities::eGaussKronrodLegendre:
419 {
420 hasEndPoints = false;
421 }
422 break;
429 {
430 hasEndPoints = true;
431 }
432 break;
433 default:
434 ASSERTL0(false, "FwdTransBndConstrained cannot be used "
435 "with this point type");
436 }
437
438 fill(outarray.data(), outarray.data() + m_ncoeffs, 0.0);
439
440 if (hasEndPoints && hasEndModes)
441 {
442
443 outarray[GetVertexMap(0)] = inarray[0];
444 outarray[GetVertexMap(1)] = inarray[m_base[0]->GetNumPoints() - 1];
445
446 if (m_ncoeffs > 2)
447 {
448 // ideally, we would like to have tmp0 to be replaced
449 // by outarray (currently MassMatrixOp does not allow
450 // aliasing)
451 Array<OneD, NekDouble> tmp0(m_ncoeffs);
452 Array<OneD, NekDouble> tmp1(m_ncoeffs);
453
454 StdRegions::StdMatrixKey stdmasskey(StdRegions::eMass,
455 DetShapeType(), *this);
456 MassMatrixOp(outarray, tmp0, stdmasskey);
457 v_IProductWRTBase(inarray, tmp1);
458
459 Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
460
461 // get Mass matrix inverse (only of interior DOF)
462 MatrixKey masskey(StdRegions::eMass, DetShapeType(), *this);
463 DNekScalMatSharedPtr matsys =
464 (m_staticCondMatrixManager[masskey])->GetBlock(1, 1);
465
466 Blas::Dgemv('N', nInteriorDofs, nInteriorDofs, matsys->Scale(),
467 &((matsys->GetOwnedMatrix())->GetPtr())[0],
468 nInteriorDofs, tmp1.data() + offset, 1, 0.0,
469 outarray.data() + offset, 1);
470 }
471 }
472 else
473 {
474 SegExp::v_FwdTrans(inarray, outarray);
475 }
476 }
477}
void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Forward transform from physical quadrature space stored in inarray and evaluate the expansion coeffic...
Definition SegExp.cpp:346
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
int GetVertexMap(const int localVertexId, bool useCoeffPacking=false)
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
static void Dgemv(const char &trans, const int &m, const int &n, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = alpha A x plus beta y where A[m x n].
Definition Blas.hpp:211
@ eGaussLegendreWithMP
1D Gauss-Legendre quadrature points with additional x=-1 and x=1 end points
Definition PointsType.h:95
@ eGaussLobattoChebyshev
1D Gauss-Lobatto-Legendre quadrature points
Definition PointsType.h:57
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
Definition PointsType.h:74
@ eGaussGaussChebyshev
1D Gauss-Gauss-Chebyshev quadrature points
Definition PointsType.h:52
@ ePolyEvenlySpaced
1D Evenly-spaced points using Lagrange polynomial
Definition PointsType.h:73
@ eGaussLobattoKronrodLegendre
1D Lobatto Kronrod quadrature points
Definition PointsType.h:72
@ eModified_B
Principle Modified Functions .
Definition BasisType.h:49
@ eOrtho_A
Principle Orthogonal Functions .
Definition BasisType.h:42
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition Vmath.hpp:220

References ASSERTL0, Nektar::StdRegions::StdExpansion::DetShapeType(), Blas::Dgemv(), Nektar::LibUtilities::eFourierEvenlySpaced, Nektar::LibUtilities::eGauss_Lagrange, Nektar::LibUtilities::eGaussGaussChebyshev, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::LibUtilities::eGaussLegendreWithMP, Nektar::LibUtilities::eGaussLobattoChebyshev, Nektar::LibUtilities::eGaussLobattoKronrodLegendre, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGLL_Lagrange, Nektar::StdRegions::eMass, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::LibUtilities::eOrtho_A, 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::StdExpansion.

Definition at line 1458 of file SegExp.cpp.

1459{
1460 DNekMatSharedPtr returnval;
1461
1462 switch (mkey.GetMatrixType())
1463 {
1470 returnval = Expansion1D::v_GenMatrix(mkey);
1471 break;
1472 default:
1473 returnval = StdSegExp::v_GenMatrix(mkey);
1474 break;
1475 }
1476
1477 return returnval;
1478}
DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey) override

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

◆ v_GetCoord()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 654 of file SegExp.cpp.

656{
657 int i;
658
659 ASSERTL1(Lcoords[0] >= -1.0 && Lcoords[0] <= 1.0,
660 "Local coordinates are not in region [-1,1]");
661
662 m_geom->FillGeom();
663 for (i = 0; i < m_geom->GetCoordim(); ++i)
664 {
665 coords[i] = m_geom->GetCoord(i, Lcoords);
666 }
667}
NekDouble GetCoord(const int i, const Array< OneD, const NekDouble > &Lcoord)
Given local collapsed coordinate Lcoord, return the value of physical coordinate in direction i.
Definition Geometry.h:558
void FillGeom()
Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.
Definition Geometry.h:460

References ASSERTL1, Nektar::SpatialDomains::Geometry::FillGeom(), Nektar::SpatialDomains::Geometry::GetCoord(), Nektar::SpatialDomains::Geometry::GetCoordim(), and Nektar::LocalRegions::Expansion::m_geom.

◆ v_GetCoords()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 669 of file SegExp.cpp.

672{
673 Expansion::v_GetCoords(coords_0, coords_1, coords_2);
674}
void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override

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

766{
767 LibUtilities::BasisKey bkey0(m_base[0]->GetBasisType(), 2,
768 m_base[0]->GetPointsKey());
769
771}

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

1109{
1110 return m_matrixManager[mkey];
1111}

References m_matrixManager.

◆ v_GetLocStaticCondMatrix()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1098 of file SegExp.cpp.

1099{
1100 return m_staticCondMatrixManager[mkey];
1101}

References m_staticCondMatrixManager.

◆ v_GetPhysNormals()

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

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 773 of file SegExp.cpp.

774{
775 NEKERROR(ErrorUtil::efatal, "Got to SegExp");
777}
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
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 725 of file SegExp.cpp.

726{
727 int nquad = m_base[0]->GetNumPoints();
728
729 ASSERTL1(vertex == 0 || vertex == 1, "Vertex value should be 0 or 1");
730
731 map = Array<OneD, int>(1);
732
733 map[0] = vertex == 0 ? 0 : nquad - 1;
734}

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

718{
719 NekDouble result;
720 v_GetVertexPhysVals(edge, inarray, result);
721 outarray[0] = result;
722}
void v_GetVertexPhysVals(const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray) override
Definition SegExp.cpp:677

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

680{
681 int nquad = m_base[0]->GetNumPoints();
682
684 {
685 switch (vertex)
686 {
687 case 0:
688 outarray = inarray[0];
689 break;
690 case 1:
691 outarray = inarray[nquad - 1];
692 break;
693 }
694 }
695 else
696 {
699
700 StdRegions::StdMatrixKey key(StdRegions::eInterpGauss, DetShapeType(),
701 *this, factors);
702
703 DNekScalMatSharedPtr mat_gauss = m_matrixManager[key];
704
705 outarray =
706 Blas::Ddot(nquad, mat_gauss->GetOwnedMatrix()->GetPtr().data(), 1,
707 &inarray[0], 1);
708 }
709}
static double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = .
Definition Blas.hpp:163

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

Definition at line 999 of file SegExp.cpp.

1002{
1003 int nquad = m_base[0]->GetNumPoints();
1004 const Array<TwoD, const NekDouble> &gmat =
1005 m_metricinfo->GetDerivFactors(GetPointsKeys());
1006 const NekDouble lambda = mkey.GetConstFactor(StdRegions::eFactorLambda);
1007
1008 Array<OneD, NekDouble> physValues(nquad);
1009 Array<OneD, NekDouble> dPhysValuesdx(nquad);
1010 Array<OneD, NekDouble> wsp(m_ncoeffs);
1011
1012 BwdTrans(inarray, physValues);
1013
1014 // mass matrix operation
1015 v_IProductWRTBase((m_base[0]->GetBdata()), physValues, wsp, 1);
1016
1017 // Laplacian matrix operation
1018 switch (m_geom->GetCoordim())
1019 {
1020 case 1:
1021 {
1022 PhysDeriv(physValues, dPhysValuesdx);
1023
1024 // multiply with the proper geometric factors
1025 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1026 {
1027 Vmath::Vmul(nquad, &gmat[0][0], 1, dPhysValuesdx.data(), 1,
1028 dPhysValuesdx.data(), 1);
1029 }
1030 else
1031 {
1032 Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.data(), 1);
1033 }
1034 }
1035 break;
1036 case 2:
1037 {
1038 Array<OneD, NekDouble> dPhysValuesdy(nquad);
1039
1040 PhysDeriv(physValues, dPhysValuesdx, dPhysValuesdy);
1041
1042 // multiply with the proper geometric factors
1043 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1044 {
1045 Vmath::Vmul(nquad, &gmat[0][0], 1, dPhysValuesdx.data(), 1,
1046 dPhysValuesdx.data(), 1);
1047 Vmath::Vvtvp(nquad, &gmat[1][0], 1, dPhysValuesdy.data(), 1,
1048 dPhysValuesdx.data(), 1, dPhysValuesdx.data(), 1);
1049 }
1050 else
1051 {
1052 Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.data(), 1);
1053 Blas::Daxpy(nquad, gmat[1][0], dPhysValuesdy.data(), 1,
1054 dPhysValuesdx.data(), 1);
1055 }
1056 }
1057 break;
1058 case 3:
1059 {
1060 Array<OneD, NekDouble> dPhysValuesdy(nquad);
1061 Array<OneD, NekDouble> dPhysValuesdz(nquad);
1062
1063 PhysDeriv(physValues, dPhysValuesdx, dPhysValuesdy, dPhysValuesdz);
1064
1065 // multiply with the proper geometric factors
1066 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1067 {
1068 Vmath::Vmul(nquad, &gmat[0][0], 1, dPhysValuesdx.data(), 1,
1069 dPhysValuesdx.data(), 1);
1070 Vmath::Vvtvp(nquad, &gmat[1][0], 1, dPhysValuesdy.data(), 1,
1071 dPhysValuesdx.data(), 1, dPhysValuesdx.data(), 1);
1072 Vmath::Vvtvp(nquad, &gmat[2][0], 1, dPhysValuesdz.data(), 1,
1073 dPhysValuesdx.data(), 1, dPhysValuesdx.data(), 1);
1074 }
1075 else
1076 {
1077 Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.data(), 1);
1078 Blas::Daxpy(nquad, gmat[1][0], dPhysValuesdy.data(), 1,
1079 dPhysValuesdx.data(), 1);
1080 Blas::Daxpy(nquad, gmat[2][0], dPhysValuesdz.data(), 1,
1081 dPhysValuesdx.data(), 1);
1082 }
1083 }
1084 break;
1085 default:
1086 ASSERTL0(false, "Wrong number of dimensions");
1087 break;
1088 }
1089
1090 v_IProductWRTBase(m_base[0]->GetDbdata(), dPhysValuesdx, outarray, 1);
1091 Blas::Daxpy(m_ncoeffs, lambda, wsp.data(), 1, outarray.data(), 1);
1092}
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Backward transformation from coefficient space to physical space.
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)
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition Blas.hpp:149
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
Definition Blas.hpp:135
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition Vmath.hpp:72
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition Vmath.hpp:366

References ASSERTL0, Nektar::StdRegions::StdExpansion::BwdTrans(), Blas::Daxpy(), Blas::Dscal(), Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::SpatialDomains::Geometry::GetCoordim(), 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::StdExpansion.

Definition at line 98 of file SegExp.cpp.

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

References Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetNumPoints(), 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::StdExpansion.

Definition at line 531 of file SegExp.cpp.

534{
535 int nquad0 = m_base[0]->GetNumPoints();
536 Array<OneD, const NekDouble> jac = m_metricinfo->GetJac(GetPointsKeys());
537 Array<OneD, NekDouble> tmp(nquad0);
538
539 // multiply inarray with Jacobian
540 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
541 {
542 Vmath::Vmul(nquad0, jac, 1, inarray, 1, tmp, 1);
543 }
544 else
545 {
546 Vmath::Smul(nquad0, jac[0], inarray, 1, tmp, 1);
547 }
548 StdSegExp::v_IProductWRTBase(base, tmp, outarray, coll_check);
549}

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 498 of file SegExp.cpp.

500{
501 v_IProductWRTBase(m_base[0]->GetBdata(), inarray, outarray, 1);
502}

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

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

◆ v_IProductWRTDerivBase()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 551 of file SegExp.cpp.

554{
555 ASSERTL1(dir < 3, "input dir is out of range");
556 ASSERTL1((dir == 2) ? m_geom->GetCoordim() == 3 : true,
557 "input dir is out of range");
558
559 int nquad = m_base[0]->GetNumPoints();
560 const Array<TwoD, const NekDouble> &gmat =
561 m_metricinfo->GetDerivFactors(GetPointsKeys());
562
563 Array<OneD, NekDouble> tmp1(nquad);
564
565 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
566 {
567 Vmath::Vmul(nquad, gmat[dir], 1, inarray, 1, tmp1, 1);
568 }
569 else
570 {
571 Vmath::Smul(nquad, gmat[dir][0], inarray, 1, tmp1, 1);
572 }
573
574 v_IProductWRTBase(m_base[0]->GetDbdata(), tmp1, outarray, 1);
575}

References ASSERTL1, Nektar::SpatialDomains::eDeformed, Nektar::SpatialDomains::Geometry::GetCoordim(), 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() [1/2]

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::StdExpansion.

Definition at line 901 of file SegExp.cpp.

905{
906 int nquad = m_base[0]->GetNumPoints();
907 const Array<TwoD, const NekDouble> &gmat =
908 m_metricinfo->GetDerivFactors(GetPointsKeys());
909
910 Array<OneD, NekDouble> physValues(nquad);
911 Array<OneD, NekDouble> dPhysValuesdx(nquad);
912
913 BwdTrans(inarray, physValues);
914
915 // Laplacian matrix operation
916 switch (m_geom->GetCoordim())
917 {
918 case 1:
919 {
920 PhysDeriv(physValues, dPhysValuesdx);
921
922 // multiply with the proper geometric factors
923 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
924 {
925 Vmath::Vmul(nquad, &gmat[0][0], 1, dPhysValuesdx.data(), 1,
926 dPhysValuesdx.data(), 1);
927 }
928 else
929 {
930 Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.data(), 1);
931 }
932 }
933 break;
934 case 2:
935 {
936 Array<OneD, NekDouble> dPhysValuesdy(nquad);
937
938 PhysDeriv(physValues, dPhysValuesdx, dPhysValuesdy);
939
940 // multiply with the proper geometric factors
941 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
942 {
943 Vmath::Vmul(nquad, &gmat[0][0], 1, dPhysValuesdx.data(), 1,
944 dPhysValuesdx.data(), 1);
945 Vmath::Vvtvp(nquad, &gmat[1][0], 1, dPhysValuesdy.data(), 1,
946 dPhysValuesdx.data(), 1, dPhysValuesdx.data(), 1);
947 }
948 else
949 {
950 Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.data(), 1);
951 Blas::Daxpy(nquad, gmat[1][0], dPhysValuesdy.data(), 1,
952 dPhysValuesdx.data(), 1);
953 }
954 }
955 break;
956 case 3:
957 {
958 Array<OneD, NekDouble> dPhysValuesdy(nquad);
959 Array<OneD, NekDouble> dPhysValuesdz(nquad);
960
961 PhysDeriv(physValues, dPhysValuesdx, dPhysValuesdy, dPhysValuesdz);
962
963 // multiply with the proper geometric factors
964 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
965 {
966 Vmath::Vmul(nquad, &gmat[0][0], 1, dPhysValuesdx.data(), 1,
967 dPhysValuesdx.data(), 1);
968 Vmath::Vvtvp(nquad, &gmat[1][0], 1, dPhysValuesdy.data(), 1,
969 dPhysValuesdx.data(), 1, dPhysValuesdx.data(), 1);
970 Vmath::Vvtvp(nquad, &gmat[2][0], 1, dPhysValuesdz.data(), 1,
971 dPhysValuesdx.data(), 1, dPhysValuesdx.data(), 1);
972 }
973 else
974 {
975 Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.data(), 1);
976 Blas::Daxpy(nquad, gmat[1][0], dPhysValuesdy.data(), 1,
977 dPhysValuesdx.data(), 1);
978 Blas::Daxpy(nquad, gmat[2][0], dPhysValuesdz.data(), 1,
979 dPhysValuesdx.data(), 1);
980 }
981 }
982 break;
983 default:
984 ASSERTL0(false, "Wrong number of dimensions");
985 break;
986 }
987
988 v_IProductWRTBase(m_base[0]->GetDbdata(), dPhysValuesdx, outarray, 1);
989}

References ASSERTL0, Nektar::StdRegions::StdExpansion::BwdTrans(), Blas::Daxpy(), Blas::Dscal(), Nektar::SpatialDomains::eDeformed, Nektar::SpatialDomains::Geometry::GetCoordim(), 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_LaplacianMatrixOp() [2/2]

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 991 of file SegExp.cpp.

995{
996 StdExpansion::LaplacianMatrixOp_MatFree(k1, k2, inarray, outarray, mkey);
997}

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

598{
599 NormVectorIProductWRTBase(Fvec[0], Fvec[1], outarray);
600}
void NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)

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

580{
581 int nq = m_base[0]->GetNumPoints();
582 Array<OneD, NekDouble> Fn(nq);
583
584 // @TODO: This routine no longer makes sense as a normal is not unique to an
585 // edge
586 const Array<OneD, const Array<OneD, NekDouble>> &normals =
587 GetLeftAdjacentElementExp()->GetTraceNormal(
589 Vmath::Vmul(nq, &Fx[0], 1, &normals[0][0], 1, &Fn[0], 1);
590 Vmath::Vvtvp(nq, &Fy[0], 1, &normals[1][0], 1, &Fn[0], 1, &Fn[0], 1);
591
592 v_IProductWRTBase(Fn, outarray);
593}
ExpansionSharedPtr GetLeftAdjacentElementExp() const
Definition Expansion.h:447
int GetLeftAdjacentElementTrace() const
Definition Expansion.h:460

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 779 of file SegExp.cpp.

780{
781 return 2;
782}

◆ v_NumDGBndryCoeffs()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 784 of file SegExp.cpp.

785{
786 return 2;
787}

◆ 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::StdExpansion.

Definition at line 145 of file SegExp.cpp.

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

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

Referenced by v_PhysDeriv_n().

◆ v_PhysDeriv() [2/2]

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

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

See also
StdRegions::StdExpansion::PhysDeriv

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 289 of file SegExp.cpp.

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

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

◆ v_PhysDeriv_n()

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

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

Parameters
inarrayfunction to derive
out_dnresult of the derivative operation

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 240 of file SegExp.cpp.

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

References ASSERTL0, Nektar::SpatialDomains::Geometry::GetCoordim(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_geom, Nektar::LocalRegions::Expansion::m_metricinfo, Nektar::NullNekDoubleArrayOfArray, v_PhysDeriv(), Vmath::Vadd(), Vmath::Vmul(), and Vmath::Zero().

◆ v_PhysDeriv_s()

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

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

Parameters
inarrayfunction to derive
out_dsresult of the derivative operation

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 201 of file SegExp.cpp.

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

References Nektar::SpatialDomains::eDeformed, Nektar::SpatialDomains::Geometry::GetCoordim(), 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_PhysEvalFirstDeriv()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 630 of file SegExp.cpp.

634{
635 Array<OneD, NekDouble> Lcoord(1);
636 ASSERTL0(m_geom, "m_geom not defined");
637 m_geom->GetLocCoords(coord, Lcoord);
638 return StdSegExp::v_PhysEvalFirstDeriv(Lcoord, inarray, firstOrderDerivs);
639}
NekDouble GetLocCoords(const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
Determine the local collapsed coordinates that correspond to a given Cartesian coordinate for this ge...
Definition Geometry.h:548

References ASSERTL0, Nektar::SpatialDomains::Geometry::GetLocCoords(), and Nektar::LocalRegions::Expansion::m_geom.

◆ v_PhysEvalFirstSecondDeriv()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 641 of file SegExp.cpp.

646{
647 Array<OneD, NekDouble> Lcoord(1);
648 ASSERTL0(m_geom, "m_geom not defined");
649 m_geom->GetLocCoords(coord, Lcoord);
650 return StdSegExp::v_PhysEvalFirstSecondDeriv(
651 Lcoord, inarray, firstOrderDerivs, secondOrderDerivs);
652}

References ASSERTL0, Nektar::SpatialDomains::Geometry::GetLocCoords(), and Nektar::LocalRegions::Expansion::m_geom.

◆ 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::StdExpansion.

Definition at line 619 of file SegExp.cpp.

621{
622 Array<OneD, NekDouble> Lcoord = Array<OneD, NekDouble>(1);
623
624 ASSERTL0(m_geom, "m_geom not defined");
625 m_geom->GetLocCoords(coord, Lcoord);
626
627 return StdExpansion1D::v_PhysEvaluate(Lcoord, physvals);
628}

References ASSERTL0, Nektar::SpatialDomains::Geometry::GetLocCoords(), 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::StdRegions::StdExpansion.

Definition at line 740 of file SegExp.cpp.

743{
744
745 if (dir == StdRegions::eBackwards)
746 {
747 if (&inarray[0] != &outarray[0])
748 {
749 Array<OneD, NekDouble> intmp(inarray);
750 ReverseCoeffsAndSign(intmp, outarray);
751 }
752 else
753 {
754 ReverseCoeffsAndSign(inarray, outarray);
755 }
756 }
757}
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:1488

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

614{
615 // Evaluate point in local (eta) coordinates.
616 return StdExpansion1D::v_PhysEvaluate(Lcoord, physvals);
617}

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