Nektar++
Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes | 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 ()
 
void AddNormTraceInt (const int dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
Expansion2DSharedPtr GetLeftAdjacentElementExp () const
 
Expansion2DSharedPtr GetRightAdjacentElementExp () const
 
int GetLeftAdjacentElementEdge () const
 
int GetRightAdjacentElementEdge () const
 
void SetAdjacentElementExp (int edge, Expansion2DSharedPtr &e)
 
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 ()
 
DNekScalMatSharedPtr GetLocMatrix (const LocalRegions::MatrixKey &mkey)
 
DNekScalMatSharedPtr GetLocMatrix (const StdRegions::MatrixType mtype, const StdRegions::ConstFactorMap &factors=StdRegions::NullConstFactorMap, const StdRegions::VarCoeffMap &varcoeffs=StdRegions::NullVarCoeffMap)
 
SpatialDomains::GeometrySharedPtr GetGeom () const
 
void Reset ()
 
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)
 
- 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 GetNedges () const
 This function returns the number of edges of the expansion domain. More...
 
int GetEdgeNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th edge. More...
 
int GetTotalEdgeIntNcoeffs () const
 
int GetEdgeNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th edge. More...
 
int DetCartesianDirOfEdge (const int edge)
 
const LibUtilities::BasisKey DetEdgeBasisKey (const int i) const
 
const LibUtilities::BasisKey DetFaceBasisKey (const int i, const int k) const
 
int GetFaceNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th face. More...
 
int GetFaceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th face. More...
 
int GetFaceIntNcoeffs (const int i) const
 
int GetTotalFaceIntNcoeffs () const
 
int GetTraceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th edge/face. More...
 
LibUtilities::PointsKey GetFacePointsKey (const int i, const int j) const
 
int NumBndryCoeffs (void) const
 
int NumDGBndryCoeffs (void) const
 
LibUtilities::BasisType GetEdgeBasisType (const int i) const
 This function returns the type of expansion basis on the i-th edge. More...
 
const LibUtilities::PointsKey GetNodalPointsKey () const
 This function returns the type of expansion Nodal point type if defined. More...
 
int GetNfaces () const
 This function returns the number of faces of the expansion domain. More...
 
int GetNtrace () const
 Returns the number of trace elements connected to this element. More...
 
LibUtilities::ShapeType DetShapeType () const
 This function returns the shape of the expansion domain. More...
 
std::shared_ptr< StdExpansionGetStdExp (void) const
 
std::shared_ptr< StdExpansionGetLinStdExp (void) const
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion ()
 
bool IsNodalNonTensorialExp ()
 
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space. More...
 
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space. More...
 
void FwdTrans_BndConstrained (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)
 
IndexMapValuesSharedPtr GetIndexMap (const IndexMapKey &ikey)
 
const Array< OneD, const NekDouble > & GetPhysNormals (void)
 
void SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
virtual void SetUpPhysNormals (const int edge)
 
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)
 
StdRegions::Orientation GetForient (int face)
 
StdRegions::Orientation GetEorient (int edge)
 
void SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
 
void SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
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 GetEdgeInteriorMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
 
void GetFaceNumModes (const int fid, const Orientation faceOrient, int &numModes0, int &numModes1)
 
void GetFaceInteriorMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
 
void GetEdgeToElementMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int P=-1)
 
void GetFaceToElementMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int nummodesA=-1, int nummodesB=-1)
 
void GetEdgePhysVals (const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Extract the physical values along edge edge from inarray into outarray following the local edge orientation and point distribution defined by defined in EdgeExp. More...
 
void GetEdgePhysVals (const int edge, const std::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetTracePhysVals (const int edge, const std::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetVertexPhysVals (const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray)
 
void GetEdgeInterpVals (const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetEdgeQFactors (const int edge, 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 GetFacePhysVals (const int face, const std::shared_ptr< StdExpansion > &FaceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=eNoOrientation)
 
void GetEdgePhysMap (const int edge, Array< OneD, int > &outarray)
 
void GetFacePhysMap (const int face, Array< OneD, int > &outarray)
 
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)
 
void AddRobinMassMatrix (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
void AddRobinEdgeContribution (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, Array< OneD, NekDouble > &coeffs)
 
NekDouble PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta. More...
 
virtual int v_GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
virtual const Array< OneD, const NekDouble > & v_GetPhysNormals (void)
 
virtual void v_SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
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)
 
virtual StdRegions::Orientation v_GetForient (int face)
 
virtual StdRegions::Orientation v_GetEorient (int edge)
 
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 NormalVectorGetEdgeNormal (const int edge) const
 
void ComputeEdgeNormal (const int edge)
 
void NegateEdgeNormal (const int edge)
 
bool EdgeNormalNegated (const int edge)
 
void ComputeFaceNormal (const int face)
 
void NegateFaceNormal (const int face)
 
bool FaceNormalNegated (const int face)
 
void ComputeVertexNormal (const int vertex)
 
void NegateVertexNormal (const int vertex)
 
bool VertexNormalNegated (const int vertex)
 
const NormalVectorGetFaceNormal (const int face) const
 
const NormalVectorGetVertexNormal (const int vertex) const
 
const NormalVectorGetSurfaceNormal (const int id) const
 
const LibUtilities::PointsKeyVector GetPointsKeys () const
 
Array< OneD, unsigned int > GetEdgeInverseBoundaryMap (int eid)
 
Array< OneD, unsigned int > GetFaceInverseBoundaryMap (int fid, StdRegions::Orientation faceOrient=eNoOrientation, int P1=-1, int P2=-1)
 
void GetInverseBoundaryMaps (Array< OneD, unsigned int > &vmap, Array< OneD, Array< OneD, unsigned int > > &emap, Array< OneD, Array< OneD, unsigned int > > &fmap)
 
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)
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion1D
 StdExpansion1D ()
 
 StdExpansion1D (int numcoeffs, const LibUtilities::BasisKey &Ba)
 
 StdExpansion1D (const StdExpansion1D &T)
 
virtual ~StdExpansion1D ()
 
void PhysTensorDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Evaluate the derivative \( d/d{\xi_1} \) at the physical quadrature points given by inarray and return in outarray. More...
 

Protected Member Functions

virtual DNekMatSharedPtr v_GenMatrix (const StdRegions::StdMatrixKey &mkey)
 
virtual void v_AddRobinMassMatrix (const int vert, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
virtual void v_AddRobinEdgeContribution (const int vert, const Array< OneD, const NekDouble > &primCoeffs, Array< OneD, NekDouble > &coeffs)
 
virtual NekDouble v_VectorFlux (const Array< OneD, Array< OneD, NekDouble > > &vec)
 
virtual void v_NegateVertexNormal (const int vertex)
 
virtual bool v_VertexNormalNegated (const int vertex)
 
- Protected Member Functions inherited from Nektar::LocalRegions::Expansion
void ComputeLaplacianMetric ()
 
void ComputeQuadratureMetric ()
 
void ComputeGmatcdotMF (const Array< TwoD, const NekDouble > &df, const Array< OneD, const NekDouble > &direction, Array< OneD, Array< OneD, NekDouble > > &dfdir)
 
virtual void v_MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_ComputeLaplacianMetric ()
 
virtual void v_GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3)
 
Array< OneD, NekDoublev_GetMF (const int dir, const int shapedim, const StdRegions::VarCoeffMap &varcoeffs)
 
Array< OneD, NekDoublev_GetMFDiv (const int dir, const StdRegions::VarCoeffMap &varcoeffs)
 
Array< OneD, NekDoublev_GetMFMag (const int dir, const StdRegions::VarCoeffMap &varcoeffs)
 
virtual DNekScalMatSharedPtr v_GetLocMatrix (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)
 
- 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...
 
IndexMapValuesSharedPtr CreateIndexMap (const IndexMapKey &ikey)
 Create an IndexMap which contains mapping information linking any specific element shape with either its boundaries, edges, faces, verteces, etc. More...
 
void BwdTrans_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase_SumFac (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GeneralMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
 
void LaplacianMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp_MatFree (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void HelmholtzMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
 
virtual NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion1D
virtual NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 

Protected Attributes

std::map< int, bool > m_negatedNormals
 
- Protected Attributes inherited from Nektar::LocalRegions::Expansion
SpatialDomains::GeometrySharedPtr m_geom
 
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
 
MetricMap m_metrics
 
- 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
 
LibUtilities::NekManager< IndexMapKey, IndexMapValues, IndexMapKey::opLessm_IndexMapManager
 
- Protected Attributes inherited from Nektar::StdRegions::StdExpansion1D
std::map< int, NormalVectorm_vertexNormals
 

Private Attributes

Expansion2DWeakPtr m_elementLeft
 
Expansion2DWeakPtr m_elementRight
 
int m_elementEdgeLeft
 
int m_elementEdgeRight
 

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.

62  : Expansion(pGeom),
64  {
65  m_elementEdgeLeft = -1;
66  m_elementEdgeRight = -1;
67  }
Expansion(SpatialDomains::GeometrySharedPtr pGeom)
Definition: Expansion.cpp:47

References m_elementEdgeLeft, and m_elementEdgeRight.

◆ ~Expansion1D()

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

Definition at line 69 of file Expansion1D.h.

69 {}

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 303 of file Expansion1D.cpp.

305  {
306  int i,n;
307  int nbndry = NumBndryCoeffs();
308  int nquad = GetNumPoints(0);
309  int ncoeffs = GetNcoeffs();
310  int coordim = GetCoordim();
311  Array<OneD, unsigned int> vmap;
312 
313  ASSERTL0(&inarray[0] != &outarray[0],"Input and output arrays use the same memory");
314 
315 
316  const Array<OneD, const NekDouble> &Basis = GetBasis(0)->GetBdata();
318 
319  GetBoundaryMap(vmap);
320 
321  // Add F = \tau <phi_i,phi_j> (note phi_i is zero if phi_j is non-zero)
322  for(i = 0; i < nbndry; ++i)
323  {
324  outarray[vmap[i]] += tau*Basis[(vmap[i]+1)*nquad-1]*Basis[(vmap[i]+1)*nquad-1]*inarray[vmap[i]];
325  outarray[vmap[i]] += tau*Basis[vmap[i]*nquad]*Basis[vmap[i]*nquad]*inarray[vmap[i]];
326  }
327 
328 
329  //===============================================================
330  // Add -\sum_i D_i^T M^{-1} G_i + E_i M^{-1} G_i =
331  // \sum_i D_i M^{-1} G_i term
332 
336  Array<OneD, NekDouble> tmpcoeff(ncoeffs,0.0);
337  DNekVec Coeffs (ncoeffs,outarray,eWrapper);
338  DNekVec Tmpcoeff(ncoeffs,tmpcoeff,eWrapper);
339 
340  for(n = 0; n < coordim; ++n)
341  {
342  // evaluate M^{-1} G
343  for(i = 0; i < ncoeffs; ++i)
344  {
345  // lower boundary (negative normal)
346  tmpcoeff[i] -= invMass(i,vmap[0])*Basis[vmap[0]*nquad]*Basis[vmap[0]*nquad]*inarray[vmap[0]];
347 
348  // upper boundary (positive normal)
349  tmpcoeff[i] += invMass(i,vmap[1])*Basis[(vmap[1]+1)*nquad-1]*Basis[(vmap[1]+1)*nquad-1]*inarray[vmap[1]];
350 
351  }
352 
353  DNekScalMat &Dmat = *GetLocMatrix(DerivType[n]);
354  Coeffs = Coeffs + Dmat*Tmpcoeff;
355  }
356  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
DNekScalMatSharedPtr GetLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:85
void GetBoundaryMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:812
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:117
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:228
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 281 of file Expansion1D.cpp.

284  {
285  boost::ignore_unused(dir);
286 
287  int k;
288  int nbndry = NumBndryCoeffs();
289  int nquad = GetNumPoints(0);
290  const Array<OneD, const NekDouble> &Basis = GetBasis(0)->GetBdata();
291  Array<OneD, unsigned int> vmap;
292 
293  GetBoundaryMap(vmap);
294 
295  // add G \lambda term (can assume G is diagonal since one
296  // of the basis is zero at boundary otherwise)
297  for(k = 0; k < nbndry; ++k)
298  {
299  outarray[vmap[k]] += (Basis[(vmap[k]+1)*nquad-1]*Basis[(vmap[k]+1)*nquad-1] - Basis[vmap[k]*nquad]*Basis[vmap[k]*nquad])*inarray[vmap[k]];
300  }
301  }

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

◆ GetGeom1D()

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

Definition at line 172 of file Expansion1D.h.

174  {
175  return std::dynamic_pointer_cast<SpatialDomains
176  ::Geometry1D>(m_geom);
177  }
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:127

Referenced by Nektar::CompressibleFlowSystem::GetElmtMinHP(), and Nektar::DiffusionLDGNS::v_InitObject().

◆ GetLeftAdjacentElementEdge()

int Nektar::LocalRegions::Expansion1D::GetLeftAdjacentElementEdge ( ) const
inline

Definition at line 143 of file Expansion1D.h.

144  {
145  return m_elementEdgeLeft;
146  }

References m_elementEdgeLeft.

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

◆ GetLeftAdjacentElementExp()

Expansion2DSharedPtr Nektar::LocalRegions::Expansion1D::GetLeftAdjacentElementExp ( ) const
inline

Definition at line 126 of file Expansion1D.h.

128  {
129  ASSERTL1(m_elementLeft.lock().get(),
130  "Left adjacent element not set.");
131  return m_elementLeft.lock();
132  }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
Expansion2DWeakPtr m_elementLeft
Definition: Expansion1D.h:119

References ASSERTL1, and m_elementLeft.

Referenced by Nektar::MultiRegions::ExpList1D::v_GetNormals(), and Nektar::LocalRegions::SegExp::v_NormVectorIProductWRTBase().

◆ GetRightAdjacentElementEdge()

int Nektar::LocalRegions::Expansion1D::GetRightAdjacentElementEdge ( ) const
inline

Definition at line 148 of file Expansion1D.h.

149  {
150  return m_elementEdgeRight;
151  }

References m_elementEdgeRight.

◆ GetRightAdjacentElementExp()

Expansion2DSharedPtr Nektar::LocalRegions::Expansion1D::GetRightAdjacentElementExp ( ) const
inline

Definition at line 134 of file Expansion1D.h.

136  {
137  ASSERTL1(m_elementLeft.lock().get(),
138  "Right adjacent element not set.");
139 
140  return m_elementRight.lock();
141  }
Expansion2DWeakPtr m_elementRight
Definition: Expansion1D.h:120

References ASSERTL1, m_elementLeft, and m_elementRight.

◆ SetAdjacentElementExp()

void Nektar::LocalRegions::Expansion1D::SetAdjacentElementExp ( int  edge,
Expansion2DSharedPtr e 
)
inline

Definition at line 153 of file Expansion1D.h.

156  {
157  if (m_elementLeft.lock().get())
158  {
159  ASSERTL1(!m_elementRight.lock().get(),
160  "Both adjacent elements already set.");
161 
162  m_elementRight = e;
163  m_elementEdgeRight = edge;
164  }
165  else
166  {
167  m_elementLeft = e;
168  m_elementEdgeLeft = edge;
169  }
170  }

References ASSERTL1, m_elementEdgeLeft, m_elementEdgeRight, m_elementLeft, and m_elementRight.

◆ v_AddRobinEdgeContribution()

void Nektar::LocalRegions::Expansion1D::v_AddRobinEdgeContribution ( const int  vert,
const Array< OneD, const NekDouble > &  primCoeffs,
Array< OneD, NekDouble > &  coeffs 
)
protectedvirtual

Given an edge and vector of element coefficients:

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 411 of file Expansion1D.cpp.

412  {
414  "Not set up for non boundary-interior expansions");
415 
416  int map = GetVertexMap(vert);
417  Vmath::Zero(GetNcoeffs(), coeffs, 1);
418  coeffs[map] = primCoeffs[0];
419  }
int GetVertexMap(const int localVertexId, bool useCoeffPacking=false)
Definition: StdExpansion.h:822
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:376

References ASSERTL1, and Vmath::Zero().

◆ v_AddRobinMassMatrix()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 358 of file Expansion1D.cpp.

359  {
360  ASSERTL0(IsBoundaryInteriorExpansion(),"Robin boundary conditions are only implemented for boundary-interior expanisons");
361  ASSERTL1(inoutmat->GetRows() == inoutmat->GetColumns(),
362  "Assuming that input matrix was square");
363 
364  // Get local Element mapping for vertex point
365  int map = GetVertexMap(vert);
366 
367  // Now need to identify a map which takes the local edge
368  // mass matrix to the matrix stored in inoutmat;
369  // This can currently be deduced from the size of the matrix
370  // - if inoutmat.m_rows() == v_NCoeffs() it is a full
371  // matrix system
372  // - if inoutmat.m_rows() == v_NumBndCoeffs() it is a
373  // boundary CG system
374 
375  int rows = inoutmat->GetRows();
376 
377  if (rows == GetNcoeffs())
378  {
379  // no need to do anything
380  }
381  else if(rows == NumBndryCoeffs()) // same as NumDGBndryCoeffs()
382  {
383  int i;
384  Array<OneD,unsigned int> bmap;
385  GetBoundaryMap(bmap);
386 
387  for(i = 0; i < 2; ++i)
388  {
389  if(map == bmap[i])
390  {
391  map = i;
392  break;
393  }
394  }
395  ASSERTL1(i != 2,"Did not find number in map");
396  }
397 
398  // assumes end points have unit magnitude
399  (*inoutmat)(map,map) += primCoeffs[0];
400 
401  }

References ASSERTL0, and ASSERTL1.

◆ v_GenMatrix()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 61 of file Expansion1D.cpp.

62  {
63  DNekMatSharedPtr returnval;
64 
65  switch(mkey.GetMatrixType())
66  {
68  {
70  "HybridDGHelmholtz matrix not set up "
71  "for non boundary-interior expansions");
72  int i;
73  NekDouble lambdaval = mkey.GetConstFactor(StdRegions::eFactorLambda);
74  NekDouble tau = mkey.GetConstFactor(StdRegions::eFactorTau);
75  int ncoeffs = GetNcoeffs();
76 
77  int coordim = GetCoordim();
78 
83  DNekMat LocMat(ncoeffs,ncoeffs);
84 
85  returnval = MemoryManager<DNekMat>::AllocateSharedPtr(ncoeffs,ncoeffs);
86  DNekMat &Mat = *returnval;
87 
88  Vmath::Zero(ncoeffs*ncoeffs,Mat.GetPtr(),1);
89 
90  for(i=0; i < coordim; ++i)
91  {
92  DNekScalMat &Dmat = *GetLocMatrix(DerivType[i]);
93 
94  Mat = Mat + Dmat*invMass*Transpose(Dmat);
95  }
96 
97  // Add end Mass Matrix Contribution
99  Mat = Mat + lambdaval*Mass;
100 
101  Array<OneD,unsigned int> bmap;
102  GetBoundaryMap(bmap);
103 
104  // Add tau*F_e using elemental mass matrices
105  for(i = 0; i < 2; ++i)
106  {
107  Mat(bmap[i],bmap[i]) = Mat(bmap[i],bmap[i]) + tau;
108  }
109  }
110  break;
112  {
113  int j,k;
114  int nbndry = NumDGBndryCoeffs();
115  int ncoeffs = GetNcoeffs();
117  factors[StdRegions::eFactorLambda] = mkey.GetConstFactor(StdRegions::eFactorLambda);
118  factors[StdRegions::eFactorTau] = mkey.GetConstFactor(StdRegions::eFactorTau);
119 
120  Array<OneD,NekDouble> lambda(nbndry);
121  DNekVec Lambda(nbndry,lambda,eWrapper);
122  Array<OneD,NekDouble> ulam(ncoeffs);
123  DNekVec Ulam(ncoeffs,ulam,eWrapper);
124  Array<OneD,NekDouble> f(ncoeffs);
125  DNekVec F(ncoeffs,f,eWrapper);
126 
127  // declare matrix space
128  returnval = MemoryManager<DNekMat>::AllocateSharedPtr(ncoeffs,nbndry);
129  DNekMat &Umat = *returnval;
130 
131  // Helmholtz matrix
133 
134  // for each degree of freedom of the lambda space
135  // calculate Umat entry
136  // Generate Lambda to U_lambda matrix
137  for(j = 0; j < nbndry; ++j)
138  {
139  Vmath::Zero(nbndry,&lambda[0],1);
140  Vmath::Zero(ncoeffs,&f[0],1);
141  lambda[j] = 1.0;
142 
144 
145  Ulam = invHmat*F; // generate Ulam from lambda
146 
147  // fill column of matrix
148  for(k = 0; k < ncoeffs; ++k)
149  {
150  Umat(k,j) = Ulam[k];
151  }
152  }
153  }
154  break;
158  {
159  int j = 0;
160  int k = 0;
161  int dir = 0;
162  int nbndry = NumDGBndryCoeffs();
163  int ncoeffs = GetNcoeffs();
164 
165  Array<OneD,NekDouble> lambda(nbndry);
166  DNekVec Lambda(nbndry,lambda,eWrapper);
167 
168  Array<OneD,NekDouble> ulam(ncoeffs);
169  DNekVec Ulam(ncoeffs,ulam,eWrapper);
170  Array<OneD,NekDouble> f(ncoeffs);
171  DNekVec F(ncoeffs,f,eWrapper);
173  factors[StdRegions::eFactorLambda] = mkey.GetConstFactor(StdRegions::eFactorLambda);
174  factors[StdRegions::eFactorTau] = mkey.GetConstFactor(StdRegions::eFactorTau);
175 
176  // declare matrix space
177  returnval = MemoryManager<DNekMat>::AllocateSharedPtr(ncoeffs,nbndry);
178  DNekMat &Qmat = *returnval;
179 
180  // Lambda to U matrix
182 
183  // Inverse mass matrix
185 
186  //Weak Derivative matrix
188  switch(mkey.GetMatrixType())
189  {
191  dir = 0;
193  break;
195  dir = 1;
197  break;
199  dir = 2;
201  break;
202  default:
203  ASSERTL0(false,"Direction not known");
204  break;
205  }
206 
207  // for each degree of freedom of the lambda space
208  // calculate Qmat entry
209  // Generate Lambda to Q_lambda matrix
210  for(j = 0; j < nbndry; ++j)
211  {
212  Vmath::Zero(nbndry,&lambda[0],1);
213  lambda[j] = 1.0;
214 
215  // for lambda[j] = 1 this is the solution to ulam
216  for(k = 0; k < ncoeffs; ++k)
217  {
218  Ulam[k] = lamToU(k,j);
219  }
220 
221 
222  // -D^T ulam
223  Vmath::Neg(ncoeffs,&ulam[0],1);
224  F = Transpose(*Dmat)*Ulam;
225 
226  // + \tilde{G} \lambda
227  AddNormTraceInt(dir,lambda,f);
228 
229  // multiply by inverse mass matrix
230  Ulam = invMass*F;
231 
232  // fill column of matrix (Qmat is in column major format)
233  Vmath::Vcopy(ncoeffs,&ulam[0],1,&(Qmat.GetPtr())[0]+j*ncoeffs,1);
234  }
235  }
236  break;
238  {
239  int j;
240  int nbndry = NumBndryCoeffs();
241 
243  factors[StdRegions::eFactorLambda] = mkey.GetConstFactor(StdRegions::eFactorLambda);
244  factors[StdRegions::eFactorTau] = mkey.GetConstFactor(StdRegions::eFactorTau);
245 
246  Array<OneD,unsigned int> bmap;
247  Array<OneD, NekDouble> lam(2);
248  GetBoundaryMap(bmap);
249 
250  // declare matrix space
251  returnval = MemoryManager<DNekMat>::AllocateSharedPtr(nbndry, nbndry);
252  DNekMat &BndMat = *returnval;
253 
254  // Matrix to map Lambda to U
256 
257  // Matrix to map Lambda to Q
259 
260  lam[0] = 1.0; lam[1] = 0.0;
261  for(j = 0; j < nbndry; ++j)
262  {
263  BndMat(0,j) = -LamToQ(bmap[0],j) - factors[StdRegions::eFactorTau]*(LamToU(bmap[0],j) - lam[j]);
264  }
265 
266  lam[0] = 0.0; lam[1] = 1.0;
267  for(j = 0; j < nbndry; ++j)
268  {
269  BndMat(1,j) = LamToQ(bmap[1],j) - factors[StdRegions::eFactorTau]*(LamToU(bmap[1],j) - lam[j]);
270  }
271  }
272  break;
273  default:
274  ASSERTL0(false,"This matrix type cannot be generated from this class");
275  break;
276  }
277 
278  return returnval;
279  }
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:294
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
NekMatrix< NekDouble, StandardMatrixTag > DNekMat
Definition: NekTypeDefs.hpp:51
NekMatrix< InnerMatrixType, BlockMatrixTag > Transpose(NekMatrix< InnerMatrixType, BlockMatrixTag > &rhs)
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
double NekDouble
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:399
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064

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

void Nektar::LocalRegions::Expansion1D::v_NegateVertexNormal ( const int  vertex)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 46 of file Expansion1D.cpp.

47  {
48  m_negatedNormals[vertex] = true;
49  for (int i = 0; i < GetCoordim(); ++i)
50  {
51  Vmath::Neg(m_vertexNormals[vertex][i].num_elements(),
52  m_vertexNormals[vertex][i], 1);
53  }
54  }
std::map< int, bool > m_negatedNormals
Definition: Expansion1D.h:96
std::map< int, NormalVector > m_vertexNormals

References Vmath::Neg().

◆ v_VectorFlux()

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

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 421 of file Expansion1D.cpp.

423  {
424  const Array<OneD, const Array<OneD, NekDouble> >
425  &normals = GetLeftAdjacentElementExp()->
427 
428  int nq = m_base[0]->GetNumPoints();
429  Array<OneD, NekDouble > Fn(nq);
430  Vmath::Vmul (nq, &vec[0][0], 1, &normals[0][0], 1, &Fn[0], 1);
431  Vmath::Vvtvp(nq, &vec[1][0], 1, &normals[1][0], 1, &Fn[0], 1, &Fn[0], 1);
432 
433  return Integral(Fn);
434  }
Expansion2DSharedPtr GetLeftAdjacentElementExp() const
Definition: Expansion1D.h:127
const NormalVector & GetEdgeNormal(const int edge) const
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:579
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:186
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:445

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

◆ v_VertexNormalNegated()

bool Nektar::LocalRegions::Expansion1D::v_VertexNormalNegated ( const int  vertex)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 56 of file Expansion1D.cpp.

57  {
58  return m_negatedNormals[vertex];
59  }

Member Data Documentation

◆ m_elementEdgeLeft

int Nektar::LocalRegions::Expansion1D::m_elementEdgeLeft
private

Definition at line 121 of file Expansion1D.h.

Referenced by Expansion1D(), GetLeftAdjacentElementEdge(), and SetAdjacentElementExp().

◆ m_elementEdgeRight

int Nektar::LocalRegions::Expansion1D::m_elementEdgeRight
private

Definition at line 122 of file Expansion1D.h.

Referenced by Expansion1D(), GetRightAdjacentElementEdge(), and SetAdjacentElementExp().

◆ m_elementLeft

Expansion2DWeakPtr Nektar::LocalRegions::Expansion1D::m_elementLeft
private

◆ m_elementRight

Expansion2DWeakPtr Nektar::LocalRegions::Expansion1D::m_elementRight
private

Definition at line 120 of file Expansion1D.h.

Referenced by GetRightAdjacentElementExp(), and SetAdjacentElementExp().

◆ m_negatedNormals

std::map<int, bool> Nektar::LocalRegions::Expansion1D::m_negatedNormals
protected

Definition at line 96 of file Expansion1D.h.