Nektar++
Public Member Functions | Protected Member Functions | List of all members
Nektar::LocalRegions::Expansion1D Class Reference

#include <Expansion1D.h>

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

Public Member Functions

 Expansion1D (SpatialDomains::Geometry1DSharedPtr pGeom)
 
virtual ~Expansion1D () override=default
 
void AddNormTraceInt (const int dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void AddHDGHelmholtzTraceTerms (const NekDouble tau, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
SpatialDomains::Geometry1DSharedPtr GetGeom1D () const
 
- Public Member Functions inherited from Nektar::LocalRegions::Expansion
 Expansion (SpatialDomains::GeometrySharedPtr pGeom)
 
 Expansion (const Expansion &pSrc)
 
virtual ~Expansion ()
 
void SetTraceExp (const int traceid, ExpansionSharedPtr &f)
 
ExpansionSharedPtr GetTraceExp (const int traceid)
 
DNekScalMatSharedPtr GetLocMatrix (const LocalRegions::MatrixKey &mkey)
 
void DropLocMatrix (const LocalRegions::MatrixKey &mkey)
 
DNekScalMatSharedPtr GetLocMatrix (const StdRegions::MatrixType mtype, const StdRegions::ConstFactorMap &factors=StdRegions::NullConstFactorMap, const StdRegions::VarCoeffMap &varcoeffs=StdRegions::NullVarCoeffMap)
 
SpatialDomains::GeometrySharedPtr GetGeom () const
 
void Reset ()
 
IndexMapValuesSharedPtr CreateIndexMap (const IndexMapKey &ikey)
 
DNekScalBlkMatSharedPtr CreateStaticCondMatrix (const MatrixKey &mkey)
 
const SpatialDomains::GeomFactorsSharedPtrGetMetricInfo () const
 
DNekMatSharedPtr BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
 
DNekMatSharedPtr BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd)
 
void ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType)
 
void AddEdgeNormBoundaryInt (const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
void AddEdgeNormBoundaryInt (const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
void AddFaceNormBoundaryInt (const int face, const std::shared_ptr< Expansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
void DGDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble >> &coeffs, Array< OneD, NekDouble > &outarray)
 
NekDouble VectorFlux (const Array< OneD, Array< OneD, NekDouble >> &vec)
 
void NormalTraceDerivFactors (Array< OneD, Array< OneD, NekDouble >> &factors, Array< OneD, Array< OneD, NekDouble >> &d0factors, Array< OneD, Array< OneD, NekDouble >> &d1factors)
 
IndexMapValuesSharedPtr GetIndexMap (const IndexMapKey &ikey)
 
void AlignVectorToCollapsedDir (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
 
ExpansionSharedPtr GetLeftAdjacentElementExp () const
 
ExpansionSharedPtr GetRightAdjacentElementExp () const
 
int GetLeftAdjacentElementTrace () const
 
int GetRightAdjacentElementTrace () const
 
void SetAdjacentElementExp (int traceid, ExpansionSharedPtr &e)
 
StdRegions::Orientation GetTraceOrient (int trace)
 
void SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void DivideByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Divided by the metric jacobi and quadrature weights. More...
 
void GetTraceQFactors (const int trace, Array< OneD, NekDouble > &outarray)
 Extract the metric factors to compute the contravariant fluxes along edge edge and stores them into outarray following the local edge orientation (i.e. anticlockwise convention). More...
 
void GetTracePhysVals (const int trace, const StdRegions::StdExpansionSharedPtr &TraceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=StdRegions::eNoOrientation)
 
void GetTracePhysMap (const int edge, Array< OneD, int > &outarray)
 
void ReOrientTracePhysMap (const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1)
 
const NormalVectorGetTraceNormal (const int id)
 
void ComputeTraceNormal (const int id)
 
const Array< OneD, const NekDouble > & GetPhysNormals (void)
 
void SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
void SetUpPhysNormals (const int trace)
 
void AddRobinMassMatrix (const int traceid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
void TraceNormLen (const int traceid, NekDouble &h, NekDouble &p)
 
void AddRobinTraceContribution (const int traceid, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs)
 
const Array< OneD, const NekDouble > & GetElmtBndNormDirElmtLen (const int nbnd) const
 
void StdDerivBaseOnTraceMat (Array< OneD, DNekMatSharedPtr > &DerivMat)
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion
 StdExpansion ()
 Default Constructor. More...
 
 StdExpansion (const int numcoeffs, const int numbases, const LibUtilities::BasisKey &Ba=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bb=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bc=LibUtilities::NullBasisKey)
 Constructor. More...
 
 StdExpansion (const StdExpansion &T)
 Copy Constructor. More...
 
virtual ~StdExpansion ()
 Destructor. More...
 
int GetNumBases () const
 This function returns the number of 1D bases used in the expansion. More...
 
const Array< OneD, const LibUtilities::BasisSharedPtr > & GetBase () const
 This function gets the shared point to basis. More...
 
const LibUtilities::BasisSharedPtrGetBasis (int dir) const
 This function gets the shared point to basis in the dir direction. More...
 
int GetNcoeffs (void) const
 This function returns the total number of coefficients used in the expansion. More...
 
int GetTotPoints () const
 This function returns the total number of quadrature points used in the element. More...
 
LibUtilities::BasisType GetBasisType (const int dir) const
 This function returns the type of basis used in the dir direction. More...
 
int GetBasisNumModes (const int dir) const
 This function returns the number of expansion modes in the dir direction. More...
 
int EvalBasisNumModesMax (void) const
 This function returns the maximum number of expansion modes over all local directions. More...
 
LibUtilities::PointsType GetPointsType (const int dir) const
 This function returns the type of quadrature points used in the dir direction. More...
 
int GetNumPoints (const int dir) const
 This function returns the number of quadrature points in the dir direction. More...
 
const Array< OneD, const NekDouble > & GetPoints (const int dir) const
 This function returns a pointer to the array containing the quadrature points in dir direction. More...
 
int GetNverts () const
 This function returns the number of vertices of the expansion domain. More...
 
int GetTraceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th trace. More...
 
int GetTraceIntNcoeffs (const int i) const
 
int GetTraceNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th trace. More...
 
const LibUtilities::BasisKey GetTraceBasisKey (const int i, int k=-1) const
 This function returns the basis key belonging to the i-th trace. More...
 
LibUtilities::PointsKey GetTracePointsKey (const int i, int k=-1) const
 This function returns the basis key belonging to the i-th trace. More...
 
int NumBndryCoeffs (void) const
 
int NumDGBndryCoeffs (void) const
 
const LibUtilities::PointsKey GetNodalPointsKey () const
 This function returns the type of expansion Nodal point type if defined. More...
 
int GetNtraces () const
 Returns the number of trace elements connected to this element. More...
 
LibUtilities::ShapeType DetShapeType () const
 This function returns the shape of the expansion domain. More...
 
std::shared_ptr< StdExpansionGetStdExp () const
 
std::shared_ptr< StdExpansionGetLinStdExp (void) const
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion () const
 
bool IsNodalNonTensorialExp ()
 
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space. More...
 
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space. More...
 
void FwdTransBndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble Integral (const Array< OneD, const NekDouble > &inarray)
 This function integrates the specified function over the domain. More...
 
void FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 This function fills the array outarray with the mode-th mode of the expansion. More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 this function calculates the inner product of a given function f with the different modes of the expansion More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check)
 
void IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
int GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void SetElmtId (const int id)
 Set the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2=NullNekDouble1DArray, Array< OneD, NekDouble > &coords_3=NullNekDouble1DArray)
 this function returns the physical coordinates of the quadrature points of the expansion More...
 
void GetCoord (const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
 given the coordinates of a point of the element in the local collapsed coordinate system, this function calculates the physical coordinates of the point More...
 
DNekMatSharedPtr GetStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr GetStdStaticCondMatrix (const StdMatrixKey &mkey)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble >> &Fvec, Array< OneD, NekDouble > &outarray)
 
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
void DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
int CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
NekDouble StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
int GetCoordim ()
 
void GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
void GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
int GetVertexMap (const int localVertexId, bool useCoeffPacking=false)
 
void GetTraceToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
 
void GetTraceCoeffMap (const unsigned int traceid, Array< OneD, unsigned int > &maparray)
 
void GetElmtTraceToTraceMap (const unsigned int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
 
void GetTraceInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eForwards)
 
void GetTraceNumModes (const int tid, int &numModes0, int &numModes1, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2)
 
void MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
DNekMatSharedPtr CreateGeneralMatrix (const StdMatrixKey &mkey)
 this function generates the mass matrix \(\mathbf{M}[i][j] = \int \phi_i(\mathbf{x}) \phi_j(\mathbf{x}) d\mathbf{x}\) More...
 
void GeneralMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
void ExponentialFilter (Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
 
void LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
DNekMatSharedPtr GenMatrix (const StdMatrixKey &mkey)
 
void PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void PhysDeriv_s (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds)
 
void PhysDeriv_n (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn)
 
void PhysDirectionalDeriv (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &outarray)
 
void StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
 This function evaluates the first derivative of the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs, std::array< NekDouble, 6 > &secondOrderDerivs)
 
NekDouble PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluateBasis (const Array< OneD, const NekDouble > &coords, int mode)
 This function evaluates the basis function mode mode at a point coords of the domain. More...
 
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta. More...
 
void LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi)
 Convert local collapsed coordinates eta into local cartesian coordinate xi. More...
 
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble >> &Fvec, Array< OneD, NekDouble > &outarray)
 
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual void v_DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
NekDouble Linf (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( L_\infty\) error \( |\epsilon|_\infty = \max |u - u_{exact}|\) where \( u_{exact}\) is given by the array sol. More...
 
NekDouble L2 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( L_2\) error, \( | \epsilon |_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 dx \right]^{1/2} d\xi_1 \) where \( u_{exact}\) is given by the array sol. More...
 
NekDouble H1 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( H^1\) error, \( | \epsilon |^1_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 + \nabla(u - u_{exact})\cdot\nabla(u - u_{exact})\cdot dx \right]^{1/2} d\xi_1 \) where \( u_{exact}\) is given by the array sol. More...
 
const LibUtilities::PointsKeyVector GetPointsKeys () const
 
DNekMatSharedPtr BuildInverseTransformationMatrix (const DNekScalMatSharedPtr &m_transformationmatrix)
 
void PhysInterpToSimplexEquiSpaced (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int npset=-1)
 This function performs an interpolation from the physical space points provided at input into an array of equispaced points which are not the collapsed coordinate. So for a tetrahedron you will only get a tetrahedral number of values. More...
 
void GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true)
 This function provides the connectivity of local simplices (triangles or tets) to connect the equispaced data points provided by PhysInterpToSimplexEquiSpaced. More...
 
void EquiSpacedToCoeffs (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs a projection/interpolation from the equispaced points sometimes used in post-processing onto the coefficient space. More...
 
template<class T >
std::shared_ptr< T > as ()
 
void IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
void GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion1D
 StdExpansion1D ()
 
 StdExpansion1D (int numcoeffs, const LibUtilities::BasisKey &Ba)
 
 StdExpansion1D (const StdExpansion1D &T)
 
virtual ~StdExpansion1D () override
 
void PhysTensorDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Evaluate the derivative \( d/d{\xi_1} \) at the physical quadrature points given by inarray and return in outarray. More...
 
NekDouble BaryTensorDeriv (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
 
NekDouble BaryTensorDeriv (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs, std::array< NekDouble, 6 > &secondOrderDerivs)
 

Protected Member Functions

virtual DNekMatSharedPtr v_GenMatrix (const StdRegions::StdMatrixKey &mkey) override
 
virtual void v_AddRobinMassMatrix (const int vert, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat) override
 
virtual void v_AddRobinTraceContribution (const int vert, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs) override
 
virtual NekDouble v_VectorFlux (const Array< OneD, Array< OneD, NekDouble >> &vec) override
 
virtual void v_NormalTraceDerivFactors (Array< OneD, Array< OneD, NekDouble >> &factors, Array< OneD, Array< OneD, NekDouble >> &d0factors, Array< OneD, Array< OneD, NekDouble >> &d1factors) override
 : This method gets all of the factors which are required as part of the Gradient Jump Penalty stabilisation and involves the product of the normal and geometric factors along the element trace. More...
 
virtual void v_ReOrientTracePhysMap (const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1) override
 
virtual 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)
 
virtual 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 ()
 
virtual int v_GetCoordim () const override
 
virtual void v_GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
 
virtual DNekScalMatSharedPtr v_GetLocMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual void v_DropLocMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual DNekMatSharedPtr v_BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
 
virtual DNekMatSharedPtr v_BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd)
 
virtual void v_ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType)
 
virtual void v_AddEdgeNormBoundaryInt (const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddEdgeNormBoundaryInt (const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddFaceNormBoundaryInt (const int face, const std::shared_ptr< Expansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
virtual void v_DGDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble >> &coeffs, Array< OneD, NekDouble > &outarray)
 
virtual 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_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
virtual void v_GetTraceQFactors (const int trace, Array< OneD, NekDouble > &outarray)
 
virtual void v_GetTracePhysVals (const int trace, const StdRegions::StdExpansionSharedPtr &TraceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient)
 
virtual void v_GetTracePhysMap (const int edge, Array< OneD, int > &outarray)
 
virtual void v_ComputeTraceNormal (const int id)
 
virtual const Array< OneD, const NekDouble > & v_GetPhysNormals ()
 
virtual void v_SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
virtual void v_SetUpPhysNormals (const int id)
 
virtual void v_GenTraceExp (const int traceid, ExpansionSharedPtr &exp)
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion
DNekMatSharedPtr CreateStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr CreateStdStaticCondMatrix (const StdMatrixKey &mkey)
 Create the static condensation of a matrix when using a boundary interior decomposition. More...
 
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase_SumFac (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GeneralMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
 
void LaplacianMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp_MatFree (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void HelmholtzMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
 
virtual NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
virtual void v_GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 
virtual void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals, NekDouble &deriv, NekDouble &deriv2)
 This function performs the barycentric interpolation of the polynomial stored in coord at a point physvals using barycentric interpolation weights in direction. More...
 
template<int DIR>
NekDouble BaryEvaluateBasis (const NekDouble &coord, const int &mode)
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals)
 Helper function to pass an unused value by reference into BaryEvaluate. More...
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals, NekDouble &deriv)
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion1D
virtual NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals) override
 
virtual NekDouble v_PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
 
virtual NekDouble v_PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs, std::array< NekDouble, 6 > &secondOrderDerivs) override
 

Additional Inherited Members

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

Detailed Description

Definition at line 56 of file Expansion1D.h.

Constructor & Destructor Documentation

◆ Expansion1D()

Nektar::LocalRegions::Expansion1D::Expansion1D ( SpatialDomains::Geometry1DSharedPtr  pGeom)
inline

Definition at line 60 of file Expansion1D.h.

61  : Expansion(pGeom), StdExpansion1D()
62  {
63  }
Expansion(SpatialDomains::GeometrySharedPtr pGeom)
Definition: Expansion.cpp:47

◆ ~Expansion1D()

virtual Nektar::LocalRegions::Expansion1D::~Expansion1D ( )
overridevirtualdefault

Member Function Documentation

◆ AddHDGHelmholtzTraceTerms()

void Nektar::LocalRegions::Expansion1D::AddHDGHelmholtzTraceTerms ( const NekDouble  tau,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)

Definition at line 315 of file Expansion1D.cpp.

318 {
319  int i, n;
320  int nbndry = NumBndryCoeffs();
321  int nquad = GetNumPoints(0);
322  int ncoeffs = GetNcoeffs();
323  int coordim = GetCoordim();
324  Array<OneD, unsigned int> vmap;
325 
326  ASSERTL0(&inarray[0] != &outarray[0],
327  "Input and output arrays use the same memory");
328 
329  const Array<OneD, const NekDouble> &Basis = GetBasis(0)->GetBdata();
331 
332  GetBoundaryMap(vmap);
333 
334  // Add F = \tau <phi_i,phi_j> (note phi_i is zero if phi_j is non-zero)
335  for (i = 0; i < nbndry; ++i)
336  {
337  outarray[vmap[i]] += tau * Basis[(vmap[i] + 1) * nquad - 1] *
338  Basis[(vmap[i] + 1) * nquad - 1] *
339  inarray[vmap[i]];
340  outarray[vmap[i]] += tau * Basis[vmap[i] * nquad] *
341  Basis[vmap[i] * nquad] * inarray[vmap[i]];
342  }
343 
344  //===============================================================
345  // Add -\sum_i D_i^T M^{-1} G_i + E_i M^{-1} G_i =
346  // \sum_i D_i M^{-1} G_i term
347 
351  Array<OneD, NekDouble> tmpcoeff(ncoeffs, 0.0);
352  DNekVec Coeffs(ncoeffs, outarray, eWrapper);
353  DNekVec Tmpcoeff(ncoeffs, tmpcoeff, eWrapper);
354 
355  for (n = 0; n < coordim; ++n)
356  {
357  // evaluate M^{-1} G
358  for (i = 0; i < ncoeffs; ++i)
359  {
360  // lower boundary (negative normal)
361  tmpcoeff[i] -= invMass(i, vmap[0]) * Basis[vmap[0] * nquad] *
362  Basis[vmap[0] * nquad] * inarray[vmap[0]];
363 
364  // upper boundary (positive normal)
365  tmpcoeff[i] += invMass(i, vmap[1]) *
366  Basis[(vmap[1] + 1) * nquad - 1] *
367  Basis[(vmap[1] + 1) * nquad - 1] * inarray[vmap[1]];
368  }
369 
370  DNekScalMat &Dmat = *GetLocMatrix(DerivType[n]);
371  Coeffs = Coeffs + Dmat * Tmpcoeff;
372  }
373 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
DNekScalMatSharedPtr GetLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:88
void GetBoundaryMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:675
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:130
const LibUtilities::BasisSharedPtr & GetBasis(int dir) const
This function gets the shared point to basis in the dir direction.
Definition: StdExpansion.h:118
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:224
NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag > DNekScalMat
NekVector< NekDouble > DNekVec
Definition: NekTypeDefs.hpp:48

References ASSERTL0, Nektar::StdRegions::eInvMass, Nektar::StdRegions::eWeakDeriv0, Nektar::StdRegions::eWeakDeriv1, Nektar::StdRegions::eWeakDeriv2, Nektar::eWrapper, and Nektar::LibUtilities::Basis::GetBdata().

◆ AddNormTraceInt()

void Nektar::LocalRegions::Expansion1D::AddNormTraceInt ( const int  dir,
Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)

Definition at line 290 of file Expansion1D.cpp.

293 {
294  boost::ignore_unused(dir);
295 
296  int k;
297  int nbndry = NumBndryCoeffs();
298  int nquad = GetNumPoints(0);
299  const Array<OneD, const NekDouble> &Basis = GetBasis(0)->GetBdata();
300  Array<OneD, unsigned int> vmap;
301 
302  GetBoundaryMap(vmap);
303 
304  // add G \lambda term (can assume G is diagonal since one
305  // of the basis is zero at boundary otherwise)
306  for (k = 0; k < nbndry; ++k)
307  {
308  outarray[vmap[k]] += (Basis[(vmap[k] + 1) * nquad - 1] *
309  Basis[(vmap[k] + 1) * nquad - 1] -
310  Basis[vmap[k] * nquad] * Basis[vmap[k] * nquad]) *
311  inarray[vmap[k]];
312  }
313 }

References Nektar::LibUtilities::Basis::GetBdata().

◆ GetGeom1D()

SpatialDomains::Geometry1DSharedPtr Nektar::LocalRegions::Expansion1D::GetGeom1D ( ) const
inline

Definition at line 108 of file Expansion1D.h.

109 {
110  return std::dynamic_pointer_cast<SpatialDomains ::Geometry1D>(m_geom);
111 }
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:275

References Nektar::LocalRegions::Expansion::m_geom.

Referenced by Nektar::VariableConverter::SetElmtMinHP(), and Nektar::DiffusionLDGNS::v_InitObject().

◆ v_AddRobinMassMatrix()

void Nektar::LocalRegions::Expansion1D::v_AddRobinMassMatrix ( const int  vert,
const Array< OneD, const NekDouble > &  primCoeffs,
DNekMatSharedPtr inoutmat 
)
overrideprotectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 375 of file Expansion1D.cpp.

378 {
380  "Robin boundary conditions are only implemented for "
381  "boundary-interior expanisons");
382  ASSERTL1(inoutmat->GetRows() == inoutmat->GetColumns(),
383  "Assuming that input matrix was square");
384 
385  // Get local Element mapping for vertex point
386  int map = GetVertexMap(vert);
387 
388  // Now need to identify a map which takes the local edge
389  // mass matrix to the matrix stored in inoutmat;
390  // This can currently be deduced from the size of the matrix
391  // - if inoutmat.m_rows() == v_NCoeffs() it is a full
392  // matrix system
393  // - if inoutmat.m_rows() == v_NumBndCoeffs() it is a
394  // boundary CG system
395 
396  int rows = inoutmat->GetRows();
397 
398  if (rows == GetNcoeffs())
399  {
400  // no need to do anything
401  }
402  else if (rows == NumBndryCoeffs()) // same as NumDGBndryCoeffs()
403  {
404  int i;
405  Array<OneD, unsigned int> bmap;
406  GetBoundaryMap(bmap);
407 
408  for (i = 0; i < 2; ++i)
409  {
410  if (map == bmap[i])
411  {
412  map = i;
413  break;
414  }
415  }
416  ASSERTL1(i != 2, "Did not find number in map");
417  }
418 
419  // assumes end points have unit magnitude
420  (*inoutmat)(map, map) += primCoeffs[0];
421 }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
int GetVertexMap(const int localVertexId, bool useCoeffPacking=false)
Definition: StdExpansion.h:685

References ASSERTL0, and ASSERTL1.

◆ v_AddRobinTraceContribution()

void Nektar::LocalRegions::Expansion1D::v_AddRobinTraceContribution ( const int  vert,
const Array< OneD, const NekDouble > &  primCoeffs,
const Array< OneD, NekDouble > &  incoeffs,
Array< OneD, NekDouble > &  coeffs 
)
overrideprotectedvirtual

Given an edge and vector of element coefficients:

  • maps those elemental coefficients corresponding to the trace into an vector.
  • update the element coefficients
  • multiplies the edge vector by the edge mass matrix
  • maps the edge coefficients back onto the elemental coefficients

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 431 of file Expansion1D.cpp.

434 {
436  "Not set up for non boundary-interior expansions");
437 
438  int map = GetVertexMap(vert);
439  coeffs[map] += primCoeffs[0] * incoeffs[map];
440 }

References ASSERTL1.

◆ v_GenMatrix()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 47 of file Expansion1D.cpp.

48 {
49  DNekMatSharedPtr returnval;
50 
51  switch (mkey.GetMatrixType())
52  {
54  {
56  "HybridDGHelmholtz matrix not set up "
57  "for non boundary-interior expansions");
58  int i;
59  NekDouble lambdaval =
60  mkey.GetConstFactor(StdRegions::eFactorLambda);
61  NekDouble tau = mkey.GetConstFactor(StdRegions::eFactorTau);
62  int ncoeffs = GetNcoeffs();
63 
64  int coordim = GetCoordim();
65 
70  DNekMat LocMat(ncoeffs, ncoeffs);
71 
72  returnval =
74  DNekMat &Mat = *returnval;
75 
76  Vmath::Zero(ncoeffs * ncoeffs, Mat.GetPtr(), 1);
77 
78  for (i = 0; i < coordim; ++i)
79  {
80  DNekScalMat &Dmat = *GetLocMatrix(DerivType[i]);
81 
82  Mat = Mat + Dmat * invMass * Transpose(Dmat);
83  }
84 
85  // Add end Mass Matrix Contribution
87  Mat = Mat + lambdaval * Mass;
88 
89  Array<OneD, unsigned int> bmap;
90  GetBoundaryMap(bmap);
91 
92  // Add tau*F_e using elemental mass matrices
93  for (i = 0; i < 2; ++i)
94  {
95  Mat(bmap[i], bmap[i]) = Mat(bmap[i], bmap[i]) + tau;
96  }
97  }
98  break;
100  {
101  int j, k;
102  int nbndry = NumDGBndryCoeffs();
103  int ncoeffs = GetNcoeffs();
105  factors[StdRegions::eFactorLambda] =
106  mkey.GetConstFactor(StdRegions::eFactorLambda);
107  factors[StdRegions::eFactorTau] =
108  mkey.GetConstFactor(StdRegions::eFactorTau);
109 
110  Array<OneD, NekDouble> lambda(nbndry);
111  DNekVec Lambda(nbndry, lambda, eWrapper);
112  Array<OneD, NekDouble> ulam(ncoeffs);
113  DNekVec Ulam(ncoeffs, ulam, eWrapper);
114  Array<OneD, NekDouble> f(ncoeffs);
115  DNekVec F(ncoeffs, f, eWrapper);
116 
117  // declare matrix space
118  returnval =
120  DNekMat &Umat = *returnval;
121 
122  // Helmholtz matrix
123  DNekScalMat &invHmat =
125 
126  // for each degree of freedom of the lambda space
127  // calculate Umat entry
128  // Generate Lambda to U_lambda matrix
129  for (j = 0; j < nbndry; ++j)
130  {
131  Vmath::Zero(nbndry, &lambda[0], 1);
132  Vmath::Zero(ncoeffs, &f[0], 1);
133  lambda[j] = 1.0;
134 
136  lambda, f);
137 
138  Ulam = invHmat * F; // generate Ulam from lambda
139 
140  // fill column of matrix
141  for (k = 0; k < ncoeffs; ++k)
142  {
143  Umat(k, j) = Ulam[k];
144  }
145  }
146  }
147  break;
151  {
152  int j = 0;
153  int k = 0;
154  int dir = 0;
155  int nbndry = NumDGBndryCoeffs();
156  int ncoeffs = GetNcoeffs();
157 
158  Array<OneD, NekDouble> lambda(nbndry);
159  DNekVec Lambda(nbndry, lambda, eWrapper);
160 
161  Array<OneD, NekDouble> ulam(ncoeffs);
162  DNekVec Ulam(ncoeffs, ulam, eWrapper);
163  Array<OneD, NekDouble> f(ncoeffs);
164  DNekVec F(ncoeffs, f, eWrapper);
166  factors[StdRegions::eFactorLambda] =
167  mkey.GetConstFactor(StdRegions::eFactorLambda);
168  factors[StdRegions::eFactorTau] =
169  mkey.GetConstFactor(StdRegions::eFactorTau);
170 
171  // declare matrix space
172  returnval =
174  DNekMat &Qmat = *returnval;
175 
176  // Lambda to U matrix
177  DNekScalMat &lamToU =
179 
180  // Inverse mass matrix
182 
183  // Weak Derivative matrix
185  switch (mkey.GetMatrixType())
186  {
188  dir = 0;
190  break;
192  dir = 1;
194  break;
196  dir = 2;
198  break;
199  default:
200  ASSERTL0(false, "Direction not known");
201  break;
202  }
203 
204  // for each degree of freedom of the lambda space
205  // calculate Qmat entry
206  // Generate Lambda to Q_lambda matrix
207  for (j = 0; j < nbndry; ++j)
208  {
209  Vmath::Zero(nbndry, &lambda[0], 1);
210  lambda[j] = 1.0;
211 
212  // for lambda[j] = 1 this is the solution to ulam
213  for (k = 0; k < ncoeffs; ++k)
214  {
215  Ulam[k] = lamToU(k, j);
216  }
217 
218  // -D^T ulam
219  Vmath::Neg(ncoeffs, &ulam[0], 1);
220  F = Transpose(*Dmat) * Ulam;
221 
222  // + \tilde{G} \lambda
223  AddNormTraceInt(dir, lambda, f);
224 
225  // multiply by inverse mass matrix
226  Ulam = invMass * F;
227 
228  // fill column of matrix (Qmat is in column major format)
229  Vmath::Vcopy(ncoeffs, &ulam[0], 1,
230  &(Qmat.GetPtr())[0] + j * ncoeffs, 1);
231  }
232  }
233  break;
235  {
236  int j;
237  int nbndry = NumBndryCoeffs();
238 
240  factors[StdRegions::eFactorLambda] =
241  mkey.GetConstFactor(StdRegions::eFactorLambda);
242  factors[StdRegions::eFactorTau] =
243  mkey.GetConstFactor(StdRegions::eFactorTau);
244 
245  Array<OneD, unsigned int> bmap;
246  Array<OneD, NekDouble> lam(2);
247  GetBoundaryMap(bmap);
248 
249  // declare matrix space
250  returnval =
252  DNekMat &BndMat = *returnval;
253 
254  // Matrix to map Lambda to U
255  DNekScalMat &LamToU =
257 
258  // Matrix to map Lambda to Q
259  DNekScalMat &LamToQ =
261 
262  lam[0] = 1.0;
263  lam[1] = 0.0;
264  for (j = 0; j < nbndry; ++j)
265  {
266  BndMat(0, j) =
267  -LamToQ(bmap[0], j) - factors[StdRegions::eFactorTau] *
268  (LamToU(bmap[0], j) - lam[j]);
269  }
270 
271  lam[0] = 0.0;
272  lam[1] = 1.0;
273  for (j = 0; j < nbndry; ++j)
274  {
275  BndMat(1, j) =
276  LamToQ(bmap[1], j) - factors[StdRegions::eFactorTau] *
277  (LamToU(bmap[1], j) - lam[j]);
278  }
279  }
280  break;
281  default:
282  ASSERTL0(false,
283  "This matrix type cannot be generated from this class");
284  break;
285  }
286 
287  return returnval;
288 }
void AddHDGHelmholtzTraceTerms(const NekDouble tau, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void AddNormTraceInt(const int dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:399
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
NekMatrix< NekDouble, StandardMatrixTag > DNekMat
Definition: NekTypeDefs.hpp:50
NekMatrix< InnerMatrixType, BlockMatrixTag > Transpose(NekMatrix< InnerMatrixType, BlockMatrixTag > &rhs)
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:518
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255

References ASSERTL0, ASSERTL1, Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eFactorTau, Nektar::StdRegions::eHybridDGHelmBndLam, Nektar::StdRegions::eHybridDGHelmholtz, Nektar::StdRegions::eHybridDGLamToQ0, Nektar::StdRegions::eHybridDGLamToQ1, Nektar::StdRegions::eHybridDGLamToQ2, Nektar::StdRegions::eHybridDGLamToU, Nektar::StdRegions::eInvHybridDGHelmholtz, Nektar::StdRegions::eInvMass, Nektar::StdRegions::eMass, Nektar::StdRegions::eWeakDeriv0, Nektar::StdRegions::eWeakDeriv1, Nektar::StdRegions::eWeakDeriv2, Nektar::eWrapper, Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdMatrixKey::GetMatrixType(), Vmath::Neg(), Nektar::Transpose(), Vmath::Vcopy(), and Vmath::Zero().

Referenced by Nektar::LocalRegions::SegExp::v_GenMatrix().

◆ v_NormalTraceDerivFactors()

void Nektar::LocalRegions::Expansion1D::v_NormalTraceDerivFactors ( Array< OneD, Array< OneD, NekDouble >> &  factors,
Array< OneD, Array< OneD, NekDouble >> &  d0factors,
Array< OneD, Array< OneD, NekDouble >> &  d1factors 
)
overrideprotectedvirtual

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

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 462 of file Expansion1D.cpp.

466 {
467  boost::ignore_unused(d0factors, d1factors); // for 2D&3D shapes
468  int nquad = GetNumPoints(0);
469  Array<TwoD, const NekDouble> gmat =
470  m_metricinfo->GetDerivFactors(GetPointsKeys());
471 
472  if (factors.size() <= 2)
473  {
474  factors = Array<OneD, Array<OneD, NekDouble>>(2);
475  factors[0] = Array<OneD, NekDouble>(1);
476  factors[1] = Array<OneD, NekDouble>(1);
477  }
478 
479  // Outwards normal
480  const Array<OneD, const Array<OneD, NekDouble>> &normal_0 =
481  GetTraceNormal(0);
482  const Array<OneD, const Array<OneD, NekDouble>> &normal_1 =
483  GetTraceNormal(1);
484 
485  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
486  {
487  factors[0][0] = gmat[0][nquad - 1] * normal_0[0][0];
488  factors[1][0] = gmat[0][0] * normal_1[0][0];
489 
490  for (int n = 1; n < normal_0.size(); ++n)
491  {
492  factors[0][0] += gmat[n][0] * normal_0[n][0];
493  factors[1][0] += gmat[n][nquad - 1] * normal_1[n][0];
494  }
495  }
496  else
497  {
498  factors[0][0] = gmat[0][0] * normal_0[0][0];
499  factors[1][0] = gmat[0][0] * normal_1[0][0];
500 
501  for (int n = 1; n < normal_0.size(); ++n)
502  {
503  factors[0][0] += gmat[n][0] * normal_0[n][0];
504  factors[1][0] += gmat[n][0] * normal_1[n][0];
505  }
506  }
507 }
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:276
const NormalVector & GetTraceNormal(const int id)
Definition: Expansion.cpp:255
const LibUtilities::PointsKeyVector GetPointsKeys() const
@ eDeformed
Geometry is curved or has non-constant factors.

References Nektar::SpatialDomains::eDeformed.

◆ v_ReOrientTracePhysMap()

void Nektar::LocalRegions::Expansion1D::v_ReOrientTracePhysMap ( const StdRegions::Orientation  orient,
Array< OneD, int > &  idmap,
const int  nq0,
const int  nq1 
)
overrideprotectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 509 of file Expansion1D.cpp.

512 {
513  boost::ignore_unused(orient, nq0, nq1);
514 
515  if (idmap.size() != 1)
516  {
517  idmap = Array<OneD, int>(1);
518  }
519 
520  idmap[0] = 0;
521 }

◆ v_TraceNormLen()

void Nektar::LocalRegions::Expansion1D::v_TraceNormLen ( const int  traceid,
NekDouble h,
NekDouble p 
)
overrideprotectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 523 of file Expansion1D.cpp.

524 {
525  boost::ignore_unused(traceid);
526  h = GetGeom()->GetVertex(1)->dist(*GetGeom()->GetVertex(0));
527  p = m_ncoeffs - 1;
528 }
SpatialDomains::GeometrySharedPtr GetGeom() const
Definition: Expansion.cpp:171

References CellMLToNektar.cellml_metadata::p.

◆ v_VectorFlux()

NekDouble Nektar::LocalRegions::Expansion1D::v_VectorFlux ( const Array< OneD, Array< OneD, NekDouble >> &  vec)
overrideprotectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 442 of file Expansion1D.cpp.

444 {
445  const Array<OneD, const Array<OneD, NekDouble>> &normals =
446  GetLeftAdjacentElementExp()->GetTraceNormal(
448 
449  int nq = m_base[0]->GetNumPoints();
450  Array<OneD, NekDouble> Fn(nq);
451  Vmath::Vmul(nq, &vec[0][0], 1, &normals[0][0], 1, &Fn[0], 1);
452  Vmath::Vvtvp(nq, &vec[1][0], 1, &normals[1][0], 1, &Fn[0], 1, &Fn[0], 1);
453 
454  return Integral(Fn);
455 }
ExpansionSharedPtr GetLeftAdjacentElementExp() const
Definition: Expansion.h:443
int GetLeftAdjacentElementTrace() const
Definition: Expansion.h:456
Array< OneD, LibUtilities::BasisSharedPtr > m_base
NekDouble Integral(const Array< OneD, const NekDouble > &inarray)
This function integrates the specified function over the domain.
Definition: StdExpansion.h:479
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:209
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:574

References Vmath::Vmul(), and Vmath::Vvtvp().