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

#include <Expansion2D.h>

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

Public Member Functions

 Expansion2D (SpatialDomains::Geometry2DSharedPtr pGeom)
 
virtual ~Expansion2D () override=default
 
DNekScalMatSharedPtr CreateMatrix (const MatrixKey &mkey)
 
void SetTraceToGeomOrientation (Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, NekDouble > &inout)
 
Array< OneD, unsigned int > GetTraceInverseBoundaryMap (int eid)
 
void AddNormTraceInt (const int dir, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble >> &edgeCoeffs, Array< OneD, NekDouble > &outarray)
 
void AddNormTraceInt (const int dir, Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, NekDouble > &outarray, const StdRegions::VarCoeffMap &varcoeffs)
 
void AddEdgeBoundaryInt (const int edge, ExpansionSharedPtr &EdgeExp, Array< OneD, NekDouble > &edgePhys, Array< OneD, NekDouble > &outarray, const StdRegions::VarCoeffMap &varcoeffs=StdRegions::NullVarCoeffMap)
 
void AddHDGHelmholtzEdgeTerms (const NekDouble tau, const int edge, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, NekDouble > &edgePhys, const StdRegions::VarCoeffMap &dirForcing, Array< OneD, NekDouble > &outarray)
 
void AddHDGHelmholtzTraceTerms (const NekDouble tau, const Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, const StdRegions::VarCoeffMap &dirForcing, Array< OneD, NekDouble > &outarray)
 
SpatialDomains::Geometry2DSharedPtr GetGeom2D () const
 
void ReOrientEdgePhysMap (const int nvert, const StdRegions::Orientation orient, const int nq0, Array< OneD, int > &idmap)
 
virtual DNekMatSharedPtr v_GenMatrix (const StdRegions::StdMatrixKey &mkey) override
 
virtual void v_GenTraceExp (const int traceid, ExpansionSharedPtr &exp) override
 
- 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::StdExpansion2D
 StdExpansion2D ()
 
 StdExpansion2D (int numcoeffs, const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb)
 
 StdExpansion2D (const StdExpansion2D &T)
 
virtual ~StdExpansion2D () override
 
void PhysTensorDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d0, Array< OneD, NekDouble > &outarray_d1)
 Calculate the 2D derivative in the local tensor/collapsed coordinate at the physical points. More...
 
NekDouble Integral (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &w0, const Array< OneD, const NekDouble > &w1)
 
NekDouble BaryTensorDeriv (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
 
void BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
 
void IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
 

Protected Member Functions

virtual void v_DGDeriv (const int dir, const Array< OneD, const NekDouble > &incoeffs, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble >> &edgeCoeffs, Array< OneD, NekDouble > &out_d) override
 
virtual void v_AddEdgeNormBoundaryInt (const int edge, const ExpansionSharedPtr &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray) override
 
virtual void v_AddEdgeNormBoundaryInt (const int edge, const ExpansionSharedPtr &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray) override
 
virtual void v_AddRobinMassMatrix (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat) override
 
virtual void v_AddRobinTraceContribution (const int traceid, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs) override
 
virtual DNekMatSharedPtr v_BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd) override
 
virtual void v_ReOrientTracePhysMap (const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1) override
 
virtual void v_SetUpPhysNormals (const int edge) override
 
virtual NekDouble v_VectorFlux (const Array< OneD, Array< OneD, NekDouble >> &vec) 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 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_NormalTraceDerivFactors (Array< OneD, Array< OneD, NekDouble >> &factors, Array< OneD, Array< OneD, NekDouble >> &d0factors, Array< OneD, Array< OneD, NekDouble >> &d1factors)
 
virtual void v_AlignVectorToCollapsedDir (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
 
virtual StdRegions::Orientation v_GetTraceOrient (int trace)
 
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)
 
- 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_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::StdExpansion2D
virtual NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals) override
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
virtual NekDouble v_PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, 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 void v_BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1)=0
 
virtual void v_IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1)=0
 
virtual void v_LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
virtual void v_HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
virtual void v_GetTraceCoeffMap (const unsigned int traceid, Array< OneD, unsigned int > &maparray) override
 
virtual void v_GetElmtTraceToTraceMap (const unsigned int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient, int P, int Q) override
 Determine the mapping to re-orientate the coefficients along the element trace (assumed to align with the standard element) into the orientation of the local trace given by edgeOrient. More...
 
virtual void v_GetTraceToElementMap (const int eid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation edgeOrient=eForwards, int P=-1, int Q=-1) override
 
virtual void v_GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat) override
 

Protected Attributes

std::vector< bool > m_requireNeg
 
- 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
 

Private Member Functions

void GetPhysEdgeVarCoeffsFromElement (const int edge, ExpansionSharedPtr &EdgeExp, const Array< OneD, const NekDouble > &varcoeff, Array< OneD, NekDouble > &outarray)
 
Array< OneD, NekDoubleGetnEdgecdotMF (const int dir, const int edge, ExpansionSharedPtr &EdgeExp_e, const Array< OneD, const Array< OneD, NekDouble >> &normals, const StdRegions::VarCoeffMap &varcoeffs)
 

Detailed Description

Definition at line 59 of file Expansion2D.h.

Constructor & Destructor Documentation

◆ Expansion2D()

Nektar::LocalRegions::Expansion2D::Expansion2D ( SpatialDomains::Geometry2DSharedPtr  pGeom)

Definition at line 54 of file Expansion2D.cpp.

56 {
57 }
Expansion(SpatialDomains::GeometrySharedPtr pGeom)
Definition: Expansion.cpp:47
StdExpansion()
Default Constructor.

◆ ~Expansion2D()

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

Member Function Documentation

◆ AddEdgeBoundaryInt()

void Nektar::LocalRegions::Expansion2D::AddEdgeBoundaryInt ( const int  edge,
ExpansionSharedPtr EdgeExp,
Array< OneD, NekDouble > &  edgePhys,
Array< OneD, NekDouble > &  outarray,
const StdRegions::VarCoeffMap varcoeffs = StdRegions::NullVarCoeffMap 
)
inline

For a given edge add the \tilde{F}_1j contributions

@TODO Variable coeffs

Definition at line 836 of file Expansion2D.cpp.

841 {
842  int i;
843  int order_e = EdgeExp->GetNcoeffs();
844  int nquad_e = EdgeExp->GetNumPoints(0);
845  Array<OneD, unsigned int> map;
846  Array<OneD, int> sign;
847  Array<OneD, NekDouble> coeff(order_e);
848 
849  GetTraceToElementMap(edge, map, sign, v_GetTraceOrient(edge));
850 
854  StdRegions::VarCoeffMap::const_iterator x;
855 
856  /// @TODO Variable coeffs
857  if ((x = varcoeffs.find(VarCoeff[0])) != varcoeffs.end())
858  {
859  Array<OneD, NekDouble> work(nquad_e);
860  GetPhysEdgeVarCoeffsFromElement(edge, EdgeExp, x->second.GetValue(),
861  work);
862  Vmath::Vmul(nquad_e, work, 1, edgePhys, 1, edgePhys, 1);
863  }
864 
865  EdgeExp->IProductWRTBase(edgePhys, coeff);
866 
867  // add data to out array
868  for (i = 0; i < order_e; ++i)
869  {
870  outarray[map[i]] += sign[i] * coeff[i];
871  }
872 }
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:49
void GetPhysEdgeVarCoeffsFromElement(const int edge, ExpansionSharedPtr &EdgeExp, const Array< OneD, const NekDouble > &varcoeff, Array< OneD, NekDouble > &outarray)
virtual StdRegions::Orientation v_GetTraceOrient(int trace)
Definition: Expansion.cpp:817
void GetTraceToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
Definition: StdExpansion.h:690
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

References Nektar::StdRegions::eVarCoeffD00, Nektar::StdRegions::eVarCoeffD11, Nektar::StdRegions::eVarCoeffD22, GetPhysEdgeVarCoeffsFromElement(), Nektar::StdRegions::StdExpansion::GetTraceToElementMap(), sign, Nektar::LocalRegions::Expansion::v_GetTraceOrient(), and Vmath::Vmul().

Referenced by AddNormTraceInt().

◆ AddHDGHelmholtzEdgeTerms()

void Nektar::LocalRegions::Expansion2D::AddHDGHelmholtzEdgeTerms ( const NekDouble  tau,
const int  edge,
Array< OneD, ExpansionSharedPtr > &  EdgeExp,
Array< OneD, NekDouble > &  edgePhys,
const StdRegions::VarCoeffMap dirForcing,
Array< OneD, NekDouble > &  outarray 
)
inline

@TODO: What direction to use here??

@TODO: Document this (probably not needed)

Definition at line 907 of file Expansion2D.cpp.

911 {
912  bool mmf = (varcoeffs.find(StdRegions::eVarCoeffMF1x) != varcoeffs.end());
913  int i, j, n;
914  int nquad_e = EdgeExp[edge]->GetNumPoints(0);
915  int order_e = EdgeExp[edge]->GetNcoeffs();
916  int coordim = mmf ? 2 : GetCoordim();
917  int ncoeffs = GetNcoeffs();
918 
919  Array<OneD, NekDouble> inval(nquad_e);
920  Array<OneD, NekDouble> outcoeff(order_e);
921  Array<OneD, NekDouble> tmpcoeff(ncoeffs);
922 
923  const Array<OneD, const Array<OneD, NekDouble>> &normals =
924  GetTraceNormal(edge);
925 
926  Array<OneD, unsigned int> emap;
927  Array<OneD, int> sign;
928 
929  StdRegions::Orientation edgedir = GetTraceOrient(edge);
930 
931  DNekVec Coeffs(ncoeffs, outarray, eWrapper);
932  DNekVec Tmpcoeff(ncoeffs, tmpcoeff, eWrapper);
933 
934  GetTraceToElementMap(edge, emap, sign, edgedir);
935 
939 
943 
944  StdRegions::VarCoeffMap::const_iterator x;
945  /// @TODO: What direction to use here??
946  if ((x = varcoeffs.find(VarCoeff[0])) != varcoeffs.end())
947  {
948  Array<OneD, NekDouble> work(nquad_e);
949  GetPhysEdgeVarCoeffsFromElement(edge, EdgeExp[edge],
950  x->second.GetValue(), work);
951  Vmath::Vmul(nquad_e, work, 1, edgePhys, 1, edgePhys, 1);
952  }
953 
954  //================================================================
955  // Add F = \tau <phi_i,in_phys>
956  // Fill edge and take inner product
957  EdgeExp[edge]->IProductWRTBase(edgePhys, outcoeff);
958  // add data to out array
959  for (i = 0; i < order_e; ++i)
960  {
961  outarray[emap[i]] += sign[i] * tau * outcoeff[i];
962  }
963  //================================================================
964 
965  //===============================================================
966  // Add -\sum_i D_i^T M^{-1} G_i + E_i M^{-1} G_i =
967  // \sum_i D_i M^{-1} G_i term
968 
969  // Two independent direction
970  DNekScalMatSharedPtr invMass;
971  for (n = 0; n < coordim; ++n)
972  {
973  if (mmf)
974  {
976  Weight[StdRegions::eVarCoeffMass] = GetMFMag(n, varcoeffs);
977 
978  MatrixKey invMasskey(StdRegions::eInvMass, DetShapeType(), *this,
980 
981  invMass = GetLocMatrix(invMasskey);
982 
983  Array<OneD, NekDouble> ncdotMF_e =
984  GetnEdgecdotMF(n, edge, EdgeExp[edge], normals, varcoeffs);
985 
986  Vmath::Vmul(nquad_e, ncdotMF_e, 1, edgePhys, 1, inval, 1);
987  }
988  else
989  {
990  Vmath::Vmul(nquad_e, normals[n], 1, edgePhys, 1, inval, 1);
992  }
993 
994  // Multiply by variable coefficient
995  /// @TODO: Document this (probably not needed)
996  // StdRegions::VarCoeffMap::const_iterator x;
997  // if ((x = varcoeffs.find(VarCoeff[n])) !=
998  // varcoeffs.end())
999  // {
1000  // GetPhysEdgeVarCoeffsFromElement(edge,EdgeExp[edge],x->second,varcoeff_work);
1001  // Vmath::Vmul(nquad_e,varcoeff_work,1,EdgeExp[edge]->GetPhys(),1,EdgeExp[edge]->UpdatePhys(),1);
1002  // }
1003 
1004  EdgeExp[edge]->IProductWRTBase(inval, outcoeff);
1005 
1006  // M^{-1} G
1007  for (i = 0; i < ncoeffs; ++i)
1008  {
1009  tmpcoeff[i] = 0;
1010  for (j = 0; j < order_e; ++j)
1011  {
1012  tmpcoeff[i] += (*invMass)(i, emap[j]) * sign[j] * outcoeff[j];
1013  }
1014  }
1015 
1016  if (mmf)
1017  {
1018  StdRegions::VarCoeffMap VarCoeffDirDeriv;
1019  VarCoeffDirDeriv[StdRegions::eVarCoeffMF] =
1020  GetMF(n, coordim, varcoeffs);
1021  VarCoeffDirDeriv[StdRegions::eVarCoeffMFDiv] =
1022  GetMFDiv(n, varcoeffs);
1023 
1024  MatrixKey Dmatkey(StdRegions::eWeakDirectionalDeriv, DetShapeType(),
1026  VarCoeffDirDeriv);
1027 
1028  DNekScalMat &Dmat = *GetLocMatrix(Dmatkey);
1029 
1030  Coeffs = Coeffs + Dmat * Tmpcoeff;
1031  }
1032  else
1033  {
1034  if (varcoeffs.find(VarCoeff[n]) != varcoeffs.end())
1035  {
1036  MatrixKey mkey(DerivType[n], DetShapeType(), *this,
1037  StdRegions::NullConstFactorMap, varcoeffs);
1038 
1039  DNekScalMat &Dmat = *GetLocMatrix(mkey);
1040  Coeffs = Coeffs + Dmat * Tmpcoeff;
1041  }
1042  else
1043  {
1044  DNekScalMat &Dmat = *GetLocMatrix(DerivType[n]);
1045  Coeffs = Coeffs + Dmat * Tmpcoeff;
1046  }
1047  }
1048  }
1049 }
Array< OneD, NekDouble > GetnEdgecdotMF(const int dir, const int edge, ExpansionSharedPtr &EdgeExp_e, const Array< OneD, const Array< OneD, NekDouble >> &normals, const StdRegions::VarCoeffMap &varcoeffs)
Array< OneD, NekDouble > GetMFMag(const int dir, const StdRegions::VarCoeffMap &varcoeffs)
Definition: Expansion.cpp:714
StdRegions::Orientation GetTraceOrient(int trace)
Definition: Expansion.h:170
DNekScalMatSharedPtr GetLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:88
Array< OneD, NekDouble > GetMFDiv(const int dir, const StdRegions::VarCoeffMap &varcoeffs)
Definition: Expansion.cpp:691
const NormalVector & GetTraceNormal(const int id)
Definition: Expansion.cpp:255
Array< OneD, NekDouble > GetMF(const int dir, const int shapedim, const StdRegions::VarCoeffMap &varcoeffs)
Definition: Expansion.cpp:638
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:130
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:373
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:400
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:343
NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag > DNekScalMat
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
NekVector< NekDouble > DNekVec
Definition: NekTypeDefs.hpp:48

References Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::StdRegions::eInvMass, Nektar::StdRegions::eVarCoeffD00, Nektar::StdRegions::eVarCoeffD11, Nektar::StdRegions::eVarCoeffD22, Nektar::StdRegions::eVarCoeffMass, Nektar::StdRegions::eVarCoeffMF, Nektar::StdRegions::eVarCoeffMF1x, Nektar::StdRegions::eVarCoeffMFDiv, Nektar::StdRegions::eWeakDeriv0, Nektar::StdRegions::eWeakDeriv1, Nektar::StdRegions::eWeakDeriv2, Nektar::StdRegions::eWeakDirectionalDeriv, Nektar::eWrapper, Nektar::StdRegions::StdExpansion::GetCoordim(), Nektar::LocalRegions::Expansion::GetLocMatrix(), Nektar::LocalRegions::Expansion::GetMF(), Nektar::LocalRegions::Expansion::GetMFDiv(), Nektar::LocalRegions::Expansion::GetMFMag(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), GetnEdgecdotMF(), GetPhysEdgeVarCoeffsFromElement(), Nektar::LocalRegions::Expansion::GetTraceNormal(), Nektar::LocalRegions::Expansion::GetTraceOrient(), Nektar::StdRegions::StdExpansion::GetTraceToElementMap(), Nektar::StdRegions::NullConstFactorMap, sign, and Vmath::Vmul().

Referenced by AddHDGHelmholtzTraceTerms().

◆ AddHDGHelmholtzTraceTerms()

void Nektar::LocalRegions::Expansion2D::AddHDGHelmholtzTraceTerms ( const NekDouble  tau,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, ExpansionSharedPtr > &  EdgeExp,
const StdRegions::VarCoeffMap dirForcing,
Array< OneD, NekDouble > &  outarray 
)
inline

Definition at line 877 of file Expansion2D.cpp.

881 {
882  ASSERTL0(&inarray[0] != &outarray[0],
883  "Input and output arrays use the same memory");
884 
885  int e, cnt, order_e, nedges = GetNtraces();
886  Array<OneD, const NekDouble> tmp;
887 
888  cnt = 0;
889 
890  for (e = 0; e < nedges; ++e)
891  {
892  order_e = EdgeExp[e]->GetNcoeffs();
893  Array<OneD, NekDouble> edgeCoeffs(order_e);
894  Array<OneD, NekDouble> edgePhys(EdgeExp[e]->GetTotPoints());
895 
896  Vmath::Vcopy(order_e, tmp = inarray + cnt, 1, edgeCoeffs, 1);
897  EdgeExp[e]->BwdTrans(edgeCoeffs, edgePhys);
898  AddHDGHelmholtzEdgeTerms(tau, e, EdgeExp, edgePhys, dirForcing,
899  outarray);
900 
901  cnt += order_e;
902  }
903 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
void AddHDGHelmholtzEdgeTerms(const NekDouble tau, const int edge, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, NekDouble > &edgePhys, const StdRegions::VarCoeffMap &dirForcing, Array< OneD, NekDouble > &outarray)
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140
int GetNtraces() const
Returns the number of trace elements connected to this element.
Definition: StdExpansion.h:357
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255

References AddHDGHelmholtzEdgeTerms(), ASSERTL0, Nektar::StdRegions::StdExpansion::GetNtraces(), Nektar::StdRegions::StdExpansion::GetTotPoints(), and Vmath::Vcopy().

Referenced by v_GenMatrix().

◆ AddNormTraceInt() [1/2]

void Nektar::LocalRegions::Expansion2D::AddNormTraceInt ( const int  dir,
Array< OneD, const NekDouble > &  inarray,
Array< OneD, ExpansionSharedPtr > &  EdgeExp,
Array< OneD, NekDouble > &  outarray,
const StdRegions::VarCoeffMap varcoeffs 
)
inline

Computes the C matrix entries due to the presence of the identity matrix in Eqn. 32.

@TODO: Document this

Definition at line 748 of file Expansion2D.cpp.

753 {
754  int i, e, cnt;
755  int order_e, nquad_e;
756  int nedges = GetNtraces();
757 
758  cnt = 0;
759  for (e = 0; e < nedges; ++e)
760  {
761  order_e = EdgeExp[e]->GetNcoeffs();
762  nquad_e = EdgeExp[e]->GetNumPoints(0);
763 
764  const Array<OneD, const Array<OneD, NekDouble>> &normals =
765  GetTraceNormal(e);
766  Array<OneD, NekDouble> edgeCoeffs(order_e);
767  Array<OneD, NekDouble> edgePhys(nquad_e);
768 
769  for (i = 0; i < order_e; ++i)
770  {
771  edgeCoeffs[i] = inarray[i + cnt];
772  }
773  cnt += order_e;
774 
775  EdgeExp[e]->BwdTrans(edgeCoeffs, edgePhys);
776 
777  // Multiply by variable coefficient
778  /// @TODO: Document this
779  // StdRegions::VarCoeffType VarCoeff[3] = {StdRegions::eVarCoeffD00,
780  // StdRegions::eVarCoeffD11,
781  // StdRegions::eVarCoeffD22};
782  // StdRegions::VarCoeffMap::const_iterator x;
783  // Array<OneD, NekDouble> varcoeff_work(nquad_e);
784 
785  // if ((x = varcoeffs.find(VarCoeff[dir])) != varcoeffs.end())
786  // {
787  // GetPhysEdgeVarCoeffsFromElement(e,EdgeExp[e],x->second,varcoeff_work);
788  // Vmath::Vmul(nquad_e,varcoeff_work,1,EdgeExp[e]->GetPhys(),1,EdgeExp[e]->UpdatePhys(),1);
789  // }
790 
791  if (varcoeffs.find(StdRegions::eVarCoeffMF1x) != varcoeffs.end())
792  {
793  // MMF case
794  Array<OneD, NekDouble> ncdotMF_e =
795  GetnEdgecdotMF(dir, e, EdgeExp[e], normals, varcoeffs);
796 
797  Vmath::Vmul(nquad_e, ncdotMF_e, 1, edgePhys, 1, edgePhys, 1);
798  }
799  else
800  {
801  Vmath::Vmul(nquad_e, normals[dir], 1, edgePhys, 1, edgePhys, 1);
802  }
803 
804  AddEdgeBoundaryInt(e, EdgeExp[e], edgePhys, outarray, varcoeffs);
805  }
806 }
void AddEdgeBoundaryInt(const int edge, ExpansionSharedPtr &EdgeExp, Array< OneD, NekDouble > &edgePhys, Array< OneD, NekDouble > &outarray, const StdRegions::VarCoeffMap &varcoeffs=StdRegions::NullVarCoeffMap)

References AddEdgeBoundaryInt(), Nektar::StdRegions::eVarCoeffMF1x, GetnEdgecdotMF(), Nektar::StdRegions::StdExpansion::GetNtraces(), Nektar::LocalRegions::Expansion::GetTraceNormal(), and Vmath::Vmul().

◆ AddNormTraceInt() [2/2]

void Nektar::LocalRegions::Expansion2D::AddNormTraceInt ( const int  dir,
Array< OneD, ExpansionSharedPtr > &  EdgeExp,
Array< OneD, Array< OneD, NekDouble >> &  edgeCoeffs,
Array< OneD, NekDouble > &  outarray 
)
inline

Definition at line 808 of file Expansion2D.cpp.

812 {
813  int e;
814  int nquad_e;
815  int nedges = GetNtraces();
816 
817  for (e = 0; e < nedges; ++e)
818  {
819  nquad_e = EdgeExp[e]->GetNumPoints(0);
820 
821  Array<OneD, NekDouble> edgePhys(nquad_e);
822  const Array<OneD, const Array<OneD, NekDouble>> &normals =
823  GetTraceNormal(e);
824 
825  EdgeExp[e]->BwdTrans(edgeCoeffs[e], edgePhys);
826 
827  Vmath::Vmul(nquad_e, normals[dir], 1, edgePhys, 1, edgePhys, 1);
828 
829  AddEdgeBoundaryInt(e, EdgeExp[e], edgePhys, outarray);
830  }
831 }

References AddEdgeBoundaryInt(), Nektar::StdRegions::StdExpansion::GetNtraces(), Nektar::LocalRegions::Expansion::GetTraceNormal(), and Vmath::Vmul().

Referenced by v_DGDeriv(), and v_GenMatrix().

◆ CreateMatrix()

DNekScalMatSharedPtr Nektar::LocalRegions::Expansion2D::CreateMatrix ( const MatrixKey mkey)

Definition at line 59 of file Expansion2D.cpp.

60 {
61  DNekScalMatSharedPtr returnval;
63 
65  "Geometric information is not set up");
66 
67  switch (mkey.GetMatrixType())
68  {
69  case StdRegions::eMass:
70  {
71  if ((m_metricinfo->GetGtype() == SpatialDomains::eDeformed) ||
72  (mkey.GetNVarCoeff()))
73  {
74  NekDouble one = 1.0;
75  DNekMatSharedPtr mat = GenMatrix(mkey);
76 
77  returnval =
79  }
80  else
81  {
82  NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
83  DNekMatSharedPtr mat = GetStdMatrix(mkey);
84 
85  returnval =
87  }
88  }
89  break;
91  {
92  MatrixKey masskey(mkey, StdRegions::eMass);
93  DNekScalMat &MassMat = *GetLocMatrix(masskey);
94 
95  // Generate a local copy of traceMat
96  MatrixKey key(mkey, StdRegions::eNormDerivOnTrace);
98 
99  ASSERTL1(mkey.ConstFactorExists(StdRegions::eFactorGJP),
100  "Need to specify eFactorGJP to construct "
101  "a HelmholtzGJP matrix");
102 
103  NekDouble factor = mkey.GetConstFactor(StdRegions::eFactorGJP);
104 
105  factor /= MassMat.Scale();
106 
107  int ntot = MassMat.GetRows() * MassMat.GetColumns();
108 
109  Vmath::Svtvp(ntot, factor, &NDTraceMat->GetPtr()[0], 1,
110  MassMat.GetRawPtr(), 1, &NDTraceMat->GetPtr()[0], 1);
111 
113  MassMat.Scale(), NDTraceMat);
114  }
115  break;
117  {
118  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
119  {
120  NekDouble one = 1.0;
121  StdRegions::StdMatrixKey masskey(StdRegions::eMass,
122  DetShapeType(), *this);
123  DNekMatSharedPtr mat = GenMatrix(masskey);
124  mat->Invert();
125 
126  returnval =
128  }
129  else
130  {
131  NekDouble fac = 1.0 / (m_metricinfo->GetJac(ptsKeys))[0];
132  DNekMatSharedPtr mat = GetStdMatrix(mkey);
133 
134  returnval =
136  }
137  }
138  break;
142  {
143  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed ||
144  mkey.GetNVarCoeff())
145  {
146  NekDouble one = 1.0;
147  DNekMatSharedPtr mat = GenMatrix(mkey);
148 
149  returnval =
151  }
152  else
153  {
154  NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
155  Array<TwoD, const NekDouble> df =
156  m_metricinfo->GetDerivFactors(ptsKeys);
157  int dir = 0;
158  if (mkey.GetMatrixType() == StdRegions::eWeakDeriv0)
159  dir = 0;
160  if (mkey.GetMatrixType() == StdRegions::eWeakDeriv1)
161  dir = 1;
162  if (mkey.GetMatrixType() == StdRegions::eWeakDeriv2)
163  dir = 2;
164 
165  MatrixKey deriv0key(StdRegions::eWeakDeriv0,
166  mkey.GetShapeType(), *this);
167  MatrixKey deriv1key(StdRegions::eWeakDeriv1,
168  mkey.GetShapeType(), *this);
169 
170  DNekMat &deriv0 = *GetStdMatrix(deriv0key);
171  DNekMat &deriv1 = *GetStdMatrix(deriv1key);
172 
173  int rows = deriv0.GetRows();
174  int cols = deriv1.GetColumns();
175 
176  DNekMatSharedPtr WeakDeriv =
178  (*WeakDeriv) =
179  df[2 * dir][0] * deriv0 + df[2 * dir + 1][0] * deriv1;
180 
182  jac, WeakDeriv);
183  }
184  }
185  break;
187  {
188  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed ||
189  mkey.GetNVarCoeff())
190  {
191  NekDouble one = 1.0;
192  DNekMatSharedPtr mat = GenMatrix(mkey);
193 
194  returnval =
196  }
197  else
198  {
199  int shapedim = 2;
200 
201  // dfdirxi = tan_{xi_x} * d \xi/dx
202  // + tan_{xi_y} * d \xi/dy
203  // + tan_{xi_z} * d \xi/dz
204  // dfdireta = tan_{eta_x} * d \eta/dx
205  // + tan_{xi_y} * d \xi/dy
206  // + tan_{xi_z} * d \xi/dz
207  NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
208  Array<TwoD, const NekDouble> df =
209  m_metricinfo->GetDerivFactors(ptsKeys);
210 
211  Array<OneD, NekDouble> direction =
212  mkey.GetVarCoeff(StdRegions::eVarCoeffMF);
213 
214  // d / dx = df[0]*deriv0 + df[1]*deriv1
215  // d / dy = df[2]*deriv0 + df[3]*deriv1
216  // d / dz = df[4]*deriv0 + df[5]*deriv1
217 
218  // dfdir[dir] = e \cdot (d/dx, d/dy, d/dz)
219  // = (e^0 * df[0] + e^1 * df[2]
220  // + e^2 * df[4]) * deriv0
221  // + (e^0 * df[1] + e^1 * df[3]
222  // + e^2 * df[5]) * deriv1
223  // dfdir[dir] = e^0 * df[2 * 0 + dir]
224  // + e^1 * df[2 * 1 + dir]
225  // + e^2 * df [ 2 * 2 + dir]
226  Array<OneD, Array<OneD, NekDouble>> dfdir(shapedim);
227  Expansion::ComputeGmatcdotMF(df, direction, dfdir);
228 
229  StdRegions::VarCoeffMap dfdirxi;
230  StdRegions::VarCoeffMap dfdireta;
231 
232  dfdirxi[StdRegions::eVarCoeffWeakDeriv] = dfdir[0];
233  dfdireta[StdRegions::eVarCoeffWeakDeriv] = dfdir[1];
234 
235  MatrixKey derivxikey(StdRegions::eWeakDeriv0,
236  mkey.GetShapeType(), *this,
238  MatrixKey derivetakey(StdRegions::eWeakDeriv1,
239  mkey.GetShapeType(), *this,
241 
242  DNekMat &derivxi = *GetStdMatrix(derivxikey);
243  DNekMat &deriveta = *GetStdMatrix(derivetakey);
244 
245  int rows = derivxi.GetRows();
246  int cols = deriveta.GetColumns();
247 
248  DNekMatSharedPtr WeakDirDeriv =
250 
251  (*WeakDirDeriv) = derivxi + deriveta;
252 
253  // Add (\nabla \cdot e^k ) Mass
254  StdRegions::VarCoeffMap DiveMass;
255  DiveMass[StdRegions::eVarCoeffMass] =
256  mkey.GetVarCoeff(StdRegions::eVarCoeffMFDiv);
257  StdRegions::StdMatrixKey stdmasskey(
258  StdRegions::eMass, mkey.GetShapeType(), *this,
260 
261  DNekMatSharedPtr DiveMassmat = GetStdMatrix(stdmasskey);
262 
263  (*WeakDirDeriv) = (*WeakDirDeriv) + (*DiveMassmat);
264 
266  jac, WeakDirDeriv);
267  }
268  break;
269  }
271  {
272  if ((m_metricinfo->GetGtype() == SpatialDomains::eDeformed) ||
273  (mkey.GetNVarCoeff() > 0) ||
274  (mkey.ConstFactorExists(StdRegions::eFactorSVVCutoffRatio)))
275  {
276  NekDouble one = 1.0;
277  DNekMatSharedPtr mat = GenMatrix(mkey);
278 
279  returnval =
281  }
282  else
283  {
284  MatrixKey lap00key(StdRegions::eLaplacian00,
285  mkey.GetShapeType(), *this);
286  MatrixKey lap01key(StdRegions::eLaplacian01,
287  mkey.GetShapeType(), *this);
288  MatrixKey lap11key(StdRegions::eLaplacian11,
289  mkey.GetShapeType(), *this);
290 
291  DNekMat &lap00 = *GetStdMatrix(lap00key);
292  DNekMat &lap01 = *GetStdMatrix(lap01key);
293  DNekMat &lap11 = *GetStdMatrix(lap11key);
294 
295  NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
296  Array<TwoD, const NekDouble> gmat =
297  m_metricinfo->GetGmat(ptsKeys);
298 
299  int rows = lap00.GetRows();
300  int cols = lap00.GetColumns();
301 
302  DNekMatSharedPtr lap =
304 
305  (*lap) = gmat[0][0] * lap00 +
306  gmat[1][0] * (lap01 + Transpose(lap01)) +
307  gmat[3][0] * lap11;
308 
309  returnval =
311  }
312  }
313  break;
315  {
316  DNekMatSharedPtr mat = GenMatrix(mkey);
317 
318  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(1.0, mat);
319  }
320  break;
322  {
323  NekDouble factor = mkey.GetConstFactor(StdRegions::eFactorLambda);
324 
325  MatrixKey masskey(mkey, StdRegions::eMass);
326  DNekScalMat &MassMat = *GetLocMatrix(masskey);
327 
328  MatrixKey lapkey(mkey, StdRegions::eLaplacian);
329  DNekScalMat &LapMat = *GetLocMatrix(lapkey);
330 
331  int rows = LapMat.GetRows();
332  int cols = LapMat.GetColumns();
333 
334  DNekMatSharedPtr helm =
336 
337  NekDouble one = 1.0;
338  (*helm) = LapMat + factor * MassMat;
339 
340  returnval =
342  }
343  break;
345  {
346  MatrixKey helmkey(mkey, StdRegions::eHelmholtz);
347  DNekScalMat &HelmMat = *GetLocMatrix(helmkey);
348 
349  // Generate a local copy of traceMat
350  MatrixKey key(mkey, StdRegions::eNormDerivOnTrace);
351  DNekMatSharedPtr NDTraceMat = Expansion2D::v_GenMatrix(key);
352 
353  ASSERTL1(mkey.ConstFactorExists(StdRegions::eFactorGJP),
354  "Need to specify eFactorGJP to construct "
355  "a HelmholtzGJP matrix");
356 
357  NekDouble factor = mkey.GetConstFactor(StdRegions::eFactorGJP);
358 
359  factor /= HelmMat.Scale();
360 
361  int ntot = HelmMat.GetRows() * HelmMat.GetColumns();
362 
363  Vmath::Svtvp(ntot, factor, &NDTraceMat->GetPtr()[0], 1,
364  HelmMat.GetRawPtr(), 1, &NDTraceMat->GetPtr()[0], 1);
365 
367  HelmMat.Scale(), NDTraceMat);
368  }
369  break;
371  {
372  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
373  {
374  NekDouble one = 1.0;
375  DNekMatSharedPtr mat = GenMatrix(mkey);
376 
377  returnval =
379  }
380  else
381  {
382  NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
383  DNekMatSharedPtr mat = GetStdMatrix(mkey);
384 
385  returnval =
387  }
388  }
389  break;
393  {
394  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
395  {
396  NekDouble one = 1.0;
397  DNekMatSharedPtr mat = GenMatrix(mkey);
398 
399  returnval =
401  }
402  else
403  {
404  NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
405 
406  const Array<TwoD, const NekDouble> &df =
407  m_metricinfo->GetDerivFactors(ptsKeys);
408  int dir = 0;
409  if (mkey.GetMatrixType() == StdRegions::eIProductWRTDerivBase0)
410  dir = 0;
411  if (mkey.GetMatrixType() == StdRegions::eIProductWRTDerivBase1)
412  dir = 1;
413  if (mkey.GetMatrixType() == StdRegions::eIProductWRTDerivBase2)
414  dir = 2;
415 
416  MatrixKey iProdDeriv0Key(StdRegions::eIProductWRTDerivBase0,
417  mkey.GetShapeType(), *this);
418  MatrixKey iProdDeriv1Key(StdRegions::eIProductWRTDerivBase1,
419  mkey.GetShapeType(), *this);
420 
421  DNekMat &stdiprod0 = *GetStdMatrix(iProdDeriv0Key);
422  DNekMat &stdiprod1 = *GetStdMatrix(iProdDeriv0Key);
423 
424  int rows = stdiprod0.GetRows();
425  int cols = stdiprod1.GetColumns();
426 
427  DNekMatSharedPtr mat =
429  (*mat) =
430  df[2 * dir][0] * stdiprod0 + df[2 * dir + 1][0] * stdiprod1;
431 
432  returnval =
434  }
435  }
436  break;
437 
439  {
440  NekDouble one = 1.0;
441 
443  *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
444 
445  DNekMatSharedPtr mat = GenMatrix(hkey);
446 
447  mat->Invert();
448 
449  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one, mat);
450  }
451  break;
453  {
455  "Matrix only setup for quad elements currently");
456  DNekMatSharedPtr m_Ix;
457  Array<OneD, NekDouble> coords(1, 0.0);
458  StdRegions::ConstFactorMap factors = mkey.GetConstFactors();
459  int edge = static_cast<int>(factors[StdRegions::eFactorGaussEdge]);
460 
461  coords[0] = (edge == 0 || edge == 3) ? -1.0 : 1.0;
462 
463  m_Ix = m_base[(edge + 1) % 2]->GetI(coords);
464 
465  returnval =
467  }
468  break;
470  {
471  NekDouble one = 1.0;
472  MatrixKey helmkey(StdRegions::eHelmholtz, mkey.GetShapeType(),
473  *this, mkey.GetConstFactors(),
474  mkey.GetVarCoeffs());
475  DNekScalBlkMatSharedPtr helmStatCond =
476  GetLocStaticCondMatrix(helmkey);
477  DNekScalMatSharedPtr A = helmStatCond->GetBlock(0, 0);
479 
481  }
482  break;
483  default:
484  {
485  NekDouble one = 1.0;
486  DNekMatSharedPtr mat = GenMatrix(mkey);
487 
488  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one, mat);
489  }
490  break;
491  }
492 
493  return returnval;
494 }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:272
virtual DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey) override
DNekMatSharedPtr BuildVertexMatrix(const DNekScalMatSharedPtr &r_bnd)
Definition: Expansion.cpp:105
void ComputeGmatcdotMF(const Array< TwoD, const NekDouble > &df, const Array< OneD, const NekDouble > &direction, Array< OneD, Array< OneD, NekDouble >> &dfdir)
Definition: Expansion.cpp:608
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:276
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:609
const LibUtilities::PointsKeyVector GetPointsKeys() const
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix(const LocalRegions::MatrixKey &mkey)
Definition: StdExpansion.h:647
DNekMatSharedPtr GenMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:843
Array< OneD, LibUtilities::BasisSharedPtr > m_base
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:250
@ eNoGeomType
No type defined.
@ eDeformed
Geometry is curved or has non-constant factors.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:399
NekMatrix< NekDouble, StandardMatrixTag > DNekMat
Definition: NekTypeDefs.hpp:50
NekMatrix< InnerMatrixType, BlockMatrixTag > Transpose(NekMatrix< InnerMatrixType, BlockMatrixTag > &rhs)
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:79
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
Definition: Vmath.cpp:622

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL1, ASSERTL2, Nektar::LocalRegions::Expansion::BuildVertexMatrix(), Nektar::LocalRegions::Expansion::ComputeGmatcdotMF(), Nektar::StdRegions::StdMatrixKey::ConstFactorExists(), Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::eFactorGaussEdge, Nektar::StdRegions::eFactorGJP, Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::eHelmholtz, Nektar::StdRegions::eHelmholtzGJP, Nektar::StdRegions::eHybridDGHelmholtz, Nektar::StdRegions::eInterpGauss, Nektar::StdRegions::eInvHybridDGHelmholtz, Nektar::StdRegions::eInvLaplacianWithUnityMean, Nektar::StdRegions::eInvMass, Nektar::StdRegions::eIProductWRTBase, Nektar::StdRegions::eIProductWRTDerivBase0, Nektar::StdRegions::eIProductWRTDerivBase1, Nektar::StdRegions::eIProductWRTDerivBase2, Nektar::StdRegions::eLaplacian, Nektar::StdRegions::eLaplacian00, Nektar::StdRegions::eLaplacian01, Nektar::StdRegions::eLaplacian11, Nektar::StdRegions::eMass, Nektar::StdRegions::eMassGJP, Nektar::SpatialDomains::eNoGeomType, Nektar::StdRegions::eNormDerivOnTrace, Nektar::StdRegions::ePreconLinearSpace, Nektar::LibUtilities::eQuadrilateral, Nektar::StdRegions::eVarCoeffMass, Nektar::StdRegions::eVarCoeffMF, Nektar::StdRegions::eVarCoeffMFDiv, Nektar::StdRegions::eVarCoeffWeakDeriv, Nektar::StdRegions::eWeakDeriv0, Nektar::StdRegions::eWeakDeriv1, Nektar::StdRegions::eWeakDeriv2, Nektar::StdRegions::eWeakDirectionalDeriv, Nektar::StdRegions::StdExpansion::GenMatrix(), Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdMatrixKey::GetConstFactors(), Nektar::LocalRegions::Expansion::GetLocMatrix(), Nektar::StdRegions::StdExpansion::GetLocStaticCondMatrix(), Nektar::StdRegions::StdMatrixKey::GetMatrixType(), Nektar::StdRegions::StdMatrixKey::GetNVarCoeff(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdMatrixKey::GetShapeType(), Nektar::StdRegions::StdExpansion::GetStdMatrix(), Nektar::StdRegions::StdMatrixKey::GetVarCoeff(), Nektar::StdRegions::StdMatrixKey::GetVarCoeffs(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_metricinfo, Nektar::StdRegions::NullConstFactorMap, Vmath::Svtvp(), Nektar::Transpose(), and v_GenMatrix().

◆ GetGeom2D()

SpatialDomains::Geometry2DSharedPtr Nektar::LocalRegions::Expansion2D::GetGeom2D ( ) const
inline

Definition at line 172 of file Expansion2D.h.

173 {
174  return std::dynamic_pointer_cast<SpatialDomains::Geometry2D>(m_geom);
175 }
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:275

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

Referenced by GetTraceInverseBoundaryMap(), Nektar::LocalRegions::QuadExp::v_FwdTransBndConstrained(), Nektar::LocalRegions::TriExp::v_FwdTransBndConstrained(), and Nektar::LocalRegions::TriExp::v_GetTraceOrient().

◆ GetnEdgecdotMF()

Array< OneD, NekDouble > Nektar::LocalRegions::Expansion2D::GetnEdgecdotMF ( const int  dir,
const int  edge,
ExpansionSharedPtr EdgeExp_e,
const Array< OneD, const Array< OneD, NekDouble >> &  normals,
const StdRegions::VarCoeffMap varcoeffs 
)
private

Definition at line 2099 of file Expansion2D.cpp.

2103 {
2104  int nquad_e = EdgeExp_e->GetNumPoints(0);
2105  int coordim = GetCoordim();
2106  int nquad0 = m_base[0]->GetNumPoints();
2107  int nquad1 = m_base[1]->GetNumPoints();
2108  int nqtot = nquad0 * nquad1;
2109 
2110  StdRegions::VarCoeffType MMFCoeffs[15] = {
2119 
2120  StdRegions::VarCoeffMap::const_iterator MFdir;
2121 
2122  Array<OneD, NekDouble> ncdotMF(nqtot, 0.0);
2123  Array<OneD, NekDouble> tmp(nqtot);
2124  Array<OneD, NekDouble> tmp_e(nquad_e);
2125  for (int k = 0; k < coordim; k++)
2126  {
2127  MFdir = varcoeffs.find(MMFCoeffs[dir * 5 + k]);
2128  tmp = MFdir->second.GetValue();
2129 
2130  GetPhysEdgeVarCoeffsFromElement(edge, EdgeExp_e, tmp, tmp_e);
2131 
2132  Vmath::Vvtvp(nquad_e, &tmp_e[0], 1, &normals[k][0], 1, &ncdotMF[0], 1,
2133  &ncdotMF[0], 1);
2134  }
2135  return ncdotMF;
2136 }
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 Nektar::StdRegions::eVarCoeffMF1Div, Nektar::StdRegions::eVarCoeffMF1Mag, Nektar::StdRegions::eVarCoeffMF1x, Nektar::StdRegions::eVarCoeffMF1y, Nektar::StdRegions::eVarCoeffMF1z, Nektar::StdRegions::eVarCoeffMF2Div, Nektar::StdRegions::eVarCoeffMF2Mag, Nektar::StdRegions::eVarCoeffMF2x, Nektar::StdRegions::eVarCoeffMF2y, Nektar::StdRegions::eVarCoeffMF2z, Nektar::StdRegions::eVarCoeffMF3Div, Nektar::StdRegions::eVarCoeffMF3Mag, Nektar::StdRegions::eVarCoeffMF3x, Nektar::StdRegions::eVarCoeffMF3y, Nektar::StdRegions::eVarCoeffMF3z, Nektar::StdRegions::StdExpansion::GetCoordim(), GetPhysEdgeVarCoeffsFromElement(), Nektar::StdRegions::StdExpansion::m_base, and Vmath::Vvtvp().

Referenced by AddHDGHelmholtzEdgeTerms(), AddNormTraceInt(), and v_GenMatrix().

◆ GetPhysEdgeVarCoeffsFromElement()

void Nektar::LocalRegions::Expansion2D::GetPhysEdgeVarCoeffsFromElement ( const int  edge,
ExpansionSharedPtr EdgeExp,
const Array< OneD, const NekDouble > &  varcoeff,
Array< OneD, NekDouble > &  outarray 
)
private

Extracts the variable coefficients along an edge

Definition at line 1054 of file Expansion2D.cpp.

1058 {
1059  Array<OneD, NekDouble> tmp(GetNcoeffs());
1060  Array<OneD, NekDouble> edgetmp(EdgeExp->GetNcoeffs());
1061 
1062  // FwdTrans varcoeffs
1063  FwdTrans(varcoeff, tmp);
1064 
1065  // Map to edge
1066  Array<OneD, unsigned int> emap;
1067  Array<OneD, int> sign;
1068  StdRegions::Orientation edgedir = GetTraceOrient(edge);
1069  GetTraceToElementMap(edge, emap, sign, edgedir);
1070 
1071  for (unsigned int i = 0; i < EdgeExp->GetNcoeffs(); ++i)
1072  {
1073  edgetmp[i] = tmp[emap[i]];
1074  }
1075 
1076  // BwdTrans
1077  EdgeExp->BwdTrans(edgetmp, outarray);
1078 }
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Forward transformation from physical space to coefficient space.

References Nektar::StdRegions::StdExpansion::FwdTrans(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), Nektar::LocalRegions::Expansion::GetTraceOrient(), Nektar::StdRegions::StdExpansion::GetTraceToElementMap(), and sign.

Referenced by AddEdgeBoundaryInt(), AddHDGHelmholtzEdgeTerms(), GetnEdgecdotMF(), and v_GenMatrix().

◆ GetTraceInverseBoundaryMap()

Array< OneD, unsigned int > Nektar::LocalRegions::Expansion2D::GetTraceInverseBoundaryMap ( int  eid)

Definition at line 2021 of file Expansion2D.cpp.

2022 {
2023  int n, j;
2024  int nEdgeCoeffs;
2025  int nBndCoeffs = NumBndryCoeffs();
2026 
2027  Array<OneD, unsigned int> bmap(nBndCoeffs);
2028  GetBoundaryMap(bmap);
2029 
2030  // Map from full system to statically condensed system (i.e reverse
2031  // GetBoundaryMap)
2032  map<int, int> invmap;
2033  for (j = 0; j < nBndCoeffs; ++j)
2034  {
2035  invmap[bmap[j]] = j;
2036  }
2037 
2038  // Number of interior edge coefficients
2039  nEdgeCoeffs = GetTraceNcoeffs(eid) - 2;
2040 
2042 
2043  Array<OneD, unsigned int> edgemaparray(nEdgeCoeffs);
2044  Array<OneD, unsigned int> maparray(nEdgeCoeffs);
2045  Array<OneD, int> signarray(nEdgeCoeffs, 1);
2046  StdRegions::Orientation eOrient = geom->GetEorient(eid);
2047 
2048  // maparray is the location of the edge within the matrix
2049  GetTraceInteriorToElementMap(eid, maparray, signarray, eOrient);
2050 
2051  for (n = 0; n < nEdgeCoeffs; ++n)
2052  {
2053  edgemaparray[n] = invmap[maparray[n]];
2054  }
2055 
2056  return edgemaparray;
2057 }
SpatialDomains::Geometry2DSharedPtr GetGeom2D() const
Definition: Expansion2D.h:172
void GetBoundaryMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:675
void GetTraceInteriorToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eForwards)
Definition: StdExpansion.h:714
int GetTraceNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th trace.
Definition: StdExpansion.h:267
std::shared_ptr< Geometry2D > Geometry2DSharedPtr
Definition: Geometry.h:65

References Nektar::StdRegions::StdExpansion::GetBoundaryMap(), GetGeom2D(), Nektar::StdRegions::StdExpansion::GetTraceInteriorToElementMap(), Nektar::StdRegions::StdExpansion::GetTraceNcoeffs(), and Nektar::StdRegions::StdExpansion::NumBndryCoeffs().

◆ ReOrientEdgePhysMap()

void Nektar::LocalRegions::Expansion2D::ReOrientEdgePhysMap ( const int  nvert,
const StdRegions::Orientation  orient,
const int  nq0,
Array< OneD, int > &  idmap 
)

◆ SetTraceToGeomOrientation()

void Nektar::LocalRegions::Expansion2D::SetTraceToGeomOrientation ( Array< OneD, ExpansionSharedPtr > &  EdgeExp,
Array< OneD, NekDouble > &  inout 
)

Definition at line 729 of file Expansion2D.cpp.

731 {
732  int i, cnt = 0;
733  int nedges = GetNtraces();
734  Array<OneD, NekDouble> e_tmp;
735 
736  for (i = 0; i < nedges; ++i)
737  {
738  EdgeExp[i]->SetCoeffsToOrientation(
739  GetTraceOrient(i), e_tmp = inout + cnt, e_tmp = inout + cnt);
740  cnt += GetTraceNcoeffs(i);
741  }
742 }

References Nektar::StdRegions::StdExpansion::GetNtraces(), Nektar::StdRegions::StdExpansion::GetTraceNcoeffs(), and Nektar::LocalRegions::Expansion::GetTraceOrient().

Referenced by v_GenMatrix().

◆ v_AddEdgeNormBoundaryInt() [1/2]

void Nektar::LocalRegions::Expansion2D::v_AddEdgeNormBoundaryInt ( const int  edge,
const ExpansionSharedPtr EdgeExp,
const Array< OneD, const NekDouble > &  Fn,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Definition at line 554 of file Expansion2D.cpp.

557 {
558  int i;
559 
560  if (m_requireNeg.size() == 0)
561  {
562  int nedges = GetNtraces();
563  m_requireNeg.resize(nedges);
564 
565  for (i = 0; i < nedges; ++i)
566  {
567  m_requireNeg[i] = false;
568 
569  ExpansionSharedPtr edgeExp = m_traceExp[i].lock();
570 
571  if (edgeExp->GetRightAdjacentElementExp())
572  {
573  if (edgeExp->GetRightAdjacentElementExp()
574  ->GetGeom()
575  ->GetGlobalID() == GetGeom()->GetGlobalID())
576  {
577  m_requireNeg[i] = true;
578  }
579  }
580  }
581  }
582 
583  IndexMapKey ikey(eEdgeToElement, DetShapeType(), GetBasisNumModes(0),
584  GetBasisNumModes(1), 0, edge, GetTraceOrient(edge));
585 
587 
588  // Order of the element
589  int order_e = map->size();
590  // Order of the trace
591  int n_coeffs = EdgeExp->GetNcoeffs();
592 
593  Array<OneD, NekDouble> edgeCoeffs(n_coeffs);
594  if (n_coeffs != order_e) // Going to orthogonal space
595  {
596  EdgeExp->FwdTrans(Fn, edgeCoeffs);
597  Expansion1DSharedPtr locExp = EdgeExp->as<Expansion1D>();
598 
599  if (m_requireNeg[edge])
600  {
601  Vmath::Neg(n_coeffs, edgeCoeffs, 1);
602  }
603 
604  Array<OneD, NekDouble> coeff(n_coeffs, 0.0);
606  ((LibUtilities::BasisType)1); // 1-->Ortho_A
607  LibUtilities::BasisKey bkey_ortho(btype,
608  EdgeExp->GetBasis(0)->GetNumModes(),
609  EdgeExp->GetBasis(0)->GetPointsKey());
610  LibUtilities::BasisKey bkey(EdgeExp->GetBasis(0)->GetBasisType(),
611  EdgeExp->GetBasis(0)->GetNumModes(),
612  EdgeExp->GetBasis(0)->GetPointsKey());
613  LibUtilities::InterpCoeff1D(bkey, edgeCoeffs, bkey_ortho, coeff);
614 
615  // Cutting high frequencies
616  for (i = order_e; i < n_coeffs; i++)
617  {
618  coeff[i] = 0.0;
619  }
620 
621  LibUtilities::InterpCoeff1D(bkey_ortho, coeff, bkey, edgeCoeffs);
622 
623  StdRegions::StdMatrixKey masskey(StdRegions::eMass,
624  LibUtilities::eSegment, *EdgeExp);
625  EdgeExp->MassMatrixOp(edgeCoeffs, edgeCoeffs, masskey);
626  }
627  else
628  {
629  EdgeExp->IProductWRTBase(Fn, edgeCoeffs);
630 
631  Expansion1DSharedPtr locExp = EdgeExp->as<Expansion1D>();
632 
633  if (m_requireNeg[edge])
634  {
635  Vmath::Neg(n_coeffs, edgeCoeffs, 1);
636  }
637  }
638 
639  // Implementation for all the basis except Gauss points
640  if (EdgeExp->GetBasis(0)->GetBasisType() != LibUtilities::eGauss_Lagrange)
641  {
642  // add data to outarray if forward edge normal is outwards
643  for (i = 0; i < order_e; ++i)
644  {
645  outarray[(*map)[i].index] += (*map)[i].sign * edgeCoeffs[i];
646  }
647  }
648  else
649  {
650  int nCoeffs0, nCoeffs1;
651  int j;
652 
654  factors[StdRegions::eFactorGaussEdge] = edge;
655  StdRegions::StdMatrixKey key(StdRegions::eGaussDG, DetShapeType(),
656  *this, factors);
657 
658  DNekMatSharedPtr mat_gauss = m_stdMatrixManager[key];
659 
660  switch (edge)
661  {
662  case 0:
663  {
664  nCoeffs1 = m_base[1]->GetNumModes();
665 
666  for (i = 0; i < order_e; ++i)
667  {
668  for (j = 0; j < nCoeffs1; j++)
669  {
670  outarray[(*map)[i].index + j * order_e] +=
671  mat_gauss->GetPtr()[j] * (*map)[i].sign *
672  edgeCoeffs[i];
673  }
674  }
675  break;
676  }
677  case 1:
678  {
679  nCoeffs0 = m_base[0]->GetNumModes();
680 
681  for (i = 0; i < order_e; ++i)
682  {
683  for (j = 0; j < nCoeffs0; j++)
684  {
685  outarray[(*map)[i].index - j] +=
686  mat_gauss->GetPtr()[order_e - 1 - j] *
687  (*map)[i].sign * edgeCoeffs[i];
688  }
689  }
690  break;
691  }
692  case 2:
693  {
694  nCoeffs1 = m_base[1]->GetNumModes();
695 
696  for (i = 0; i < order_e; ++i)
697  {
698  for (j = 0; j < nCoeffs1; j++)
699  {
700  outarray[(*map)[i].index - j * order_e] +=
701  mat_gauss->GetPtr()[order_e - 1 - j] *
702  (*map)[i].sign * edgeCoeffs[i];
703  }
704  }
705  break;
706  }
707  case 3:
708  {
709  nCoeffs0 = m_base[0]->GetNumModes();
710 
711  for (i = 0; i < order_e; ++i)
712  {
713  for (j = 0; j < nCoeffs0; j++)
714  {
715  outarray[(*map)[i].index + j] +=
716  mat_gauss->GetPtr()[j] * (*map)[i].sign *
717  edgeCoeffs[i];
718  }
719  }
720  break;
721  }
722  default:
723  ASSERTL0(false, "edge value (< 3) is out of range");
724  break;
725  }
726  }
727 }
std::vector< bool > m_requireNeg
Definition: Expansion2D.h:118
SpatialDomains::GeometrySharedPtr GetGeom() const
Definition: Expansion.cpp:171
std::map< int, ExpansionWeakPtr > m_traceExp
Definition: Expansion.h:274
IndexMapValuesSharedPtr GetIndexMap(const IndexMapKey &ikey)
Definition: Expansion.h:148
LibUtilities::NekManager< StdMatrixKey, DNekMat, StdMatrixKey::opLess > m_stdMatrixManager
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:175
void InterpCoeff1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
Definition: InterpCoeff.cpp:46
@ eGauss_Lagrange
Lagrange Polynomials using the Gauss points.
Definition: BasisType.h:59
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
std::shared_ptr< IndexMapValues > IndexMapValuesSharedPtr
Definition: IndexMapKey.h:128
std::shared_ptr< Expansion1D > Expansion1DSharedPtr
Definition: Expansion1D.h:51
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:518

References ASSERTL0, Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::LocalRegions::eEdgeToElement, Nektar::StdRegions::eFactorGaussEdge, Nektar::LibUtilities::eGauss_Lagrange, Nektar::StdRegions::eGaussDG, Nektar::StdRegions::eMass, Nektar::LibUtilities::eSegment, Nektar::StdRegions::StdExpansion::GetBasisNumModes(), Nektar::LocalRegions::Expansion::GetGeom(), Nektar::LocalRegions::Expansion::GetIndexMap(), Nektar::StdRegions::StdExpansion::GetNtraces(), Nektar::LocalRegions::Expansion::GetTraceOrient(), Nektar::LibUtilities::InterpCoeff1D(), Nektar::StdRegions::StdExpansion::m_base, m_requireNeg, Nektar::StdRegions::StdExpansion::m_stdMatrixManager, Nektar::LocalRegions::Expansion::m_traceExp, and Vmath::Neg().

◆ v_AddEdgeNormBoundaryInt() [2/2]

void Nektar::LocalRegions::Expansion2D::v_AddEdgeNormBoundaryInt ( const int  edge,
const ExpansionSharedPtr EdgeExp,
const Array< OneD, const NekDouble > &  Fx,
const Array< OneD, const NekDouble > &  Fy,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Definition at line 496 of file Expansion2D.cpp.

500 {
501  ASSERTL1(GetCoordim() == 2, "Routine only set up for two-dimensions");
502 
503  const Array<OneD, const Array<OneD, NekDouble>> normals =
504  GetTraceNormal(edge);
505 
506  if (m_requireNeg.size() == 0)
507  {
508  int nedges = GetNtraces();
509  m_requireNeg.resize(nedges);
510 
511  for (int i = 0; i < nedges; ++i)
512  {
513  m_requireNeg[i] = false;
514 
515  ExpansionSharedPtr edgeExp = m_traceExp[i].lock();
516 
517  if (edgeExp->GetRightAdjacentElementExp())
518  {
519  if (edgeExp->GetRightAdjacentElementExp()
520  ->GetGeom()
521  ->GetGlobalID() == GetGeom()->GetGlobalID())
522  {
523  m_requireNeg[i] = true;
524  }
525  }
526  }
527  }
528 
529  // We allow the case of mixed polynomial order by supporting only
530  // those modes on the edge common to both adjoining elements. This
531  // is enforced here by taking the minimum size and padding with
532  // zeros.
533  int nquad_e = min(EdgeExp->GetNumPoints(0), int(normals[0].size()));
534 
535  int nEdgePts = EdgeExp->GetTotPoints();
536  Array<OneD, NekDouble> edgePhys(nEdgePts);
537  Vmath::Vmul(nquad_e, normals[0], 1, Fx, 1, edgePhys, 1);
538  Vmath::Vvtvp(nquad_e, normals[1], 1, Fy, 1, edgePhys, 1, edgePhys, 1);
539 
540  Expansion1DSharedPtr locExp = EdgeExp->as<Expansion1D>();
541 
542  if (m_requireNeg[edge])
543  {
544  if (locExp->GetRightAdjacentElementExp()->GetGeom()->GetGlobalID() ==
545  m_geom->GetGlobalID())
546  {
547  Vmath::Neg(nquad_e, edgePhys, 1);
548  }
549  }
550 
551  AddEdgeNormBoundaryInt(edge, EdgeExp, edgePhys, outarray);
552 }
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)
Definition: Expansion.cpp:118

References Nektar::LocalRegions::Expansion::AddEdgeNormBoundaryInt(), ASSERTL1, Nektar::StdRegions::StdExpansion::GetCoordim(), Nektar::LocalRegions::Expansion::GetGeom(), Nektar::StdRegions::StdExpansion::GetNtraces(), Nektar::LocalRegions::Expansion::GetTraceNormal(), Nektar::LocalRegions::Expansion::m_geom, m_requireNeg, Nektar::LocalRegions::Expansion::m_traceExp, Vmath::Neg(), Vmath::Vmul(), and Vmath::Vvtvp().

◆ v_AddRobinMassMatrix()

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

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1817 of file Expansion2D.cpp.

1820 {
1822  "Not set up for non boundary-interior expansions");
1823  ASSERTL1(inoutmat->GetRows() == inoutmat->GetColumns(),
1824  "Assuming that input matrix was square");
1825  int i, j;
1826  int id1, id2;
1827  ExpansionSharedPtr edgeExp = m_traceExp[edge].lock();
1828  int order_e = edgeExp->GetNcoeffs();
1829 
1830  Array<OneD, unsigned int> map;
1831  Array<OneD, int> sign;
1832 
1833  StdRegions::VarCoeffMap varcoeffs;
1834  varcoeffs[StdRegions::eVarCoeffMass] = primCoeffs;
1835 
1836  LocalRegions::MatrixKey mkey(StdRegions::eMass, LibUtilities::eSegment,
1838  varcoeffs);
1839  DNekScalMat &edgemat = *edgeExp->GetLocMatrix(mkey);
1840 
1841  // Now need to identify a map which takes the local edge
1842  // mass matrix to the matrix stored in inoutmat;
1843  // This can currently be deduced from the size of the matrix
1844 
1845  // - if inoutmat.m_rows() == v_NCoeffs() it is a full
1846  // matrix system
1847 
1848  // - if inoutmat.m_rows() == v_NumBndCoeffs() it is a
1849  // boundary CG system
1850 
1851  // - if inoutmat.m_rows() == v_NumDGBndCoeffs() it is a
1852  // trace DG system
1853  int rows = inoutmat->GetRows();
1854 
1855  if (rows == GetNcoeffs())
1856  {
1857  GetTraceToElementMap(edge, map, sign, v_GetTraceOrient(edge));
1858  }
1859  else if (rows == NumBndryCoeffs())
1860  {
1861  int nbndry = NumBndryCoeffs();
1862  Array<OneD, unsigned int> bmap(nbndry);
1863 
1864  GetTraceToElementMap(edge, map, sign, v_GetTraceOrient(edge));
1865 
1866  GetBoundaryMap(bmap);
1867 
1868  for (i = 0; i < order_e; ++i)
1869  {
1870  for (j = 0; j < nbndry; ++j)
1871  {
1872  if (map[i] == bmap[j])
1873  {
1874  map[i] = j;
1875  break;
1876  }
1877  }
1878  ASSERTL1(j != nbndry, "Did not find number in map");
1879  }
1880  }
1881  else if (rows == NumDGBndryCoeffs())
1882  {
1883  // possibly this should be a separate method
1884  int cnt = 0;
1885  map = Array<OneD, unsigned int>(order_e);
1886  sign = Array<OneD, int>(order_e, 1);
1887 
1888  for (i = 0; i < edge; ++i)
1889  {
1890  cnt += GetTraceNcoeffs(i);
1891  }
1892 
1893  for (i = 0; i < order_e; ++i)
1894  {
1895  map[i] = cnt++;
1896  }
1897  // check for mapping reversal
1899  {
1900  switch (edgeExp->GetBasis(0)->GetBasisType())
1901  {
1903  reverse(map.get(), map.get() + order_e);
1904  break;
1906  reverse(map.get(), map.get() + order_e);
1907  break;
1909  {
1910  swap(map[0], map[1]);
1911  for (i = 3; i < order_e; i += 2)
1912  {
1913  sign[i] = -1;
1914  }
1915  }
1916  break;
1917  default:
1918  ASSERTL0(false,
1919  "Edge boundary type not valid for this method");
1920  }
1921  }
1922  }
1923  else
1924  {
1925  ASSERTL0(false, "Could not identify matrix type from dimension");
1926  }
1927 
1928  for (i = 0; i < order_e; ++i)
1929  {
1930  id1 = map[i];
1931  for (j = 0; j < order_e; ++j)
1932  {
1933  id2 = map[j];
1934  (*inoutmat)(id1, id2) += edgemat(i, j) * sign[i] * sign[j];
1935  }
1936  }
1937 }
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:58
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:50

References ASSERTL0, ASSERTL1, Nektar::StdRegions::eBackwards, Nektar::LibUtilities::eGauss_Lagrange, Nektar::LibUtilities::eGLL_Lagrange, Nektar::StdRegions::eMass, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eSegment, Nektar::StdRegions::eVarCoeffMass, Nektar::StdRegions::StdExpansion::GetBoundaryMap(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), Nektar::StdRegions::StdExpansion::GetTraceNcoeffs(), Nektar::LocalRegions::Expansion::GetTraceOrient(), Nektar::StdRegions::StdExpansion::GetTraceToElementMap(), Nektar::StdRegions::StdExpansion::IsBoundaryInteriorExpansion(), Nektar::LocalRegions::Expansion::m_traceExp, Nektar::StdRegions::NullConstFactorMap, Nektar::StdRegions::StdExpansion::NumBndryCoeffs(), Nektar::StdRegions::StdExpansion::NumDGBndryCoeffs(), sign, and Nektar::LocalRegions::Expansion::v_GetTraceOrient().

◆ v_AddRobinTraceContribution()

void Nektar::LocalRegions::Expansion2D::v_AddRobinTraceContribution ( const int  edgeid,
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 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::LocalRegions::Expansion.

Definition at line 1947 of file Expansion2D.cpp.

1950 {
1952  "Not set up for non boundary-interior expansions");
1953  int i;
1954  ExpansionSharedPtr edgeExp = m_traceExp[edgeid].lock();
1955  int order_e = edgeExp->GetNcoeffs();
1956 
1957  Array<OneD, unsigned int> map;
1958  Array<OneD, int> sign;
1959 
1960  StdRegions::VarCoeffMap varcoeffs;
1961  varcoeffs[StdRegions::eVarCoeffMass] = primCoeffs;
1962 
1963  LocalRegions::MatrixKey mkey(StdRegions::eMass, LibUtilities::eSegment,
1965  varcoeffs);
1966  DNekScalMat &edgemat = *edgeExp->GetLocMatrix(mkey);
1967 
1968  NekVector<NekDouble> vEdgeCoeffs(order_e);
1969 
1970  GetTraceToElementMap(edgeid, map, sign, v_GetTraceOrient(edgeid));
1971 
1972  for (i = 0; i < order_e; ++i)
1973  {
1974  vEdgeCoeffs[i] = incoeffs[map[i]] * sign[i];
1975  }
1976 
1977  vEdgeCoeffs = edgemat * vEdgeCoeffs;
1978 
1979  for (i = 0; i < order_e; ++i)
1980  {
1981  coeffs[map[i]] += vEdgeCoeffs[i] * sign[i];
1982  }
1983 }

References ASSERTL1, Nektar::StdRegions::eMass, Nektar::LibUtilities::eSegment, Nektar::StdRegions::eVarCoeffMass, Nektar::StdRegions::StdExpansion::GetTraceToElementMap(), Nektar::StdRegions::StdExpansion::IsBoundaryInteriorExpansion(), Nektar::LocalRegions::Expansion::m_traceExp, Nektar::StdRegions::NullConstFactorMap, sign, and Nektar::LocalRegions::Expansion::v_GetTraceOrient().

◆ v_BuildVertexMatrix()

DNekMatSharedPtr Nektar::LocalRegions::Expansion2D::v_BuildVertexMatrix ( const DNekScalMatSharedPtr r_bnd)
overrideprotectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1985 of file Expansion2D.cpp.

1987 {
1988  MatrixStorage storage = eFULL;
1989  DNekMatSharedPtr m_vertexmatrix;
1990 
1991  int nVerts, vid1, vid2, vMap1, vMap2;
1992  NekDouble VertexValue;
1993 
1994  nVerts = GetNverts();
1995 
1996  m_vertexmatrix =
1997  MemoryManager<DNekMat>::AllocateSharedPtr(nVerts, nVerts, 0.0, storage);
1998  DNekMat &VertexMat = (*m_vertexmatrix);
1999 
2000  for (vid1 = 0; vid1 < nVerts; ++vid1)
2001  {
2002  vMap1 = GetVertexMap(vid1);
2003 
2004  for (vid2 = 0; vid2 < nVerts; ++vid2)
2005  {
2006  vMap2 = GetVertexMap(vid2);
2007  VertexValue = (*r_bnd)(vMap1, vMap2);
2008  VertexMat.SetValue(vid1, vid2, VertexValue);
2009  }
2010  }
2011 
2012  return m_vertexmatrix;
2013 }
int GetVertexMap(const int localVertexId, bool useCoeffPacking=false)
Definition: StdExpansion.h:685
int GetNverts() const
This function returns the number of vertices of the expansion domain.
Definition: StdExpansion.h:252

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::eFULL, Nektar::StdRegions::StdExpansion::GetNverts(), and Nektar::StdRegions::StdExpansion::GetVertexMap().

◆ v_DGDeriv()

void Nektar::LocalRegions::Expansion2D::v_DGDeriv ( const int  dir,
const Array< OneD, const NekDouble > &  incoeffs,
Array< OneD, ExpansionSharedPtr > &  EdgeExp,
Array< OneD, Array< OneD, NekDouble >> &  edgeCoeffs,
Array< OneD, NekDouble > &  out_d 
)
overrideprotectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1780 of file Expansion2D.cpp.

1785 {
1789 
1790  int ncoeffs = GetNcoeffs();
1791 
1793  DNekScalMat &Dmat = *GetLocMatrix(DerivType[dir]);
1794 
1795  Array<OneD, NekDouble> coeffs = incoeffs;
1796  DNekVec Coeffs(ncoeffs, coeffs, eWrapper);
1797 
1798  Coeffs = Transpose(Dmat) * Coeffs;
1799  Vmath::Neg(ncoeffs, coeffs, 1);
1800 
1801  // Add the boundary integral including the relevant part of
1802  // the normal
1803  AddNormTraceInt(dir, EdgeExp, edgeCoeffs, coeffs);
1804 
1805  DNekVec Out_d(ncoeffs, out_d, eWrapper);
1806 
1807  Out_d = InvMass * Coeffs;
1808 }
void AddNormTraceInt(const int dir, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble >> &edgeCoeffs, Array< OneD, NekDouble > &outarray)

References AddNormTraceInt(), Nektar::StdRegions::eInvMass, Nektar::StdRegions::eWeakDeriv0, Nektar::StdRegions::eWeakDeriv1, Nektar::StdRegions::eWeakDeriv2, Nektar::eWrapper, Nektar::LocalRegions::Expansion::GetLocMatrix(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), Vmath::Neg(), and Nektar::Transpose().

◆ v_GenMatrix()

DNekMatSharedPtr Nektar::LocalRegions::Expansion2D::v_GenMatrix ( const StdRegions::StdMatrixKey mkey)
overridevirtual

Computes matrices needed for the HDG formulation. References to equations relate to the following paper: R. M. Kirby, S. J. Sherwin, B. Cockburn, To CG or to HDG: A Comparative Study, J. Sci. Comp P1-30 DOI 10.1007/s10915-011-9501-7

TODO: Add variable coeffs

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp, Nektar::LocalRegions::QuadExp, and Nektar::LocalRegions::NodalTriExp.

Definition at line 1087 of file Expansion2D.cpp.

1088 {
1089  DNekMatSharedPtr returnval;
1090 
1091  switch (mkey.GetMatrixType())
1092  {
1093  // (Z^e)^{-1} (Eqn. 33, P22)
1095  {
1097  "HybridDGHelmholtz matrix not set up "
1098  "for non boundary-interior expansions");
1099 
1100  int i, j, k;
1101  NekDouble lambdaval =
1102  mkey.GetConstFactor(StdRegions::eFactorLambda);
1103  NekDouble tau = mkey.GetConstFactor(StdRegions::eFactorTau);
1104  int ncoeffs = GetNcoeffs();
1105  int nedges = GetNtraces();
1106  int shapedim = 2;
1107  const StdRegions::VarCoeffMap &varcoeffs = mkey.GetVarCoeffs();
1108  bool mmf =
1109  (varcoeffs.find(StdRegions::eVarCoeffMF1x) != varcoeffs.end());
1110 
1111  Array<OneD, unsigned int> emap;
1112  Array<OneD, int> sign;
1114  ExpansionSharedPtr EdgeExp;
1115 
1116  int order_e, coordim = GetCoordim();
1121  DNekMat LocMat(ncoeffs, ncoeffs);
1122 
1123  returnval =
1125  DNekMat &Mat = *returnval;
1126  Vmath::Zero(ncoeffs * ncoeffs, Mat.GetPtr(), 1);
1127 
1131 
1132  StdRegions::VarCoeffMap::const_iterator x;
1133 
1134  for (i = 0; i < coordim; ++i)
1135  {
1136  if (mmf)
1137  {
1138  if (i < shapedim)
1139  {
1140  StdRegions::VarCoeffMap VarCoeffDirDeriv;
1141  VarCoeffDirDeriv[StdRegions::eVarCoeffMF] =
1142  GetMF(i, shapedim, varcoeffs);
1143  VarCoeffDirDeriv[StdRegions::eVarCoeffMFDiv] =
1144  GetMFDiv(i, varcoeffs);
1145 
1146  MatrixKey Dmatkey(StdRegions::eWeakDirectionalDeriv,
1147  DetShapeType(), *this,
1149  VarCoeffDirDeriv);
1150 
1151  DNekScalMat &Dmat = *GetLocMatrix(Dmatkey);
1152 
1153  StdRegions::VarCoeffMap Weight;
1154  Weight[StdRegions::eVarCoeffMass] =
1155  GetMFMag(i, mkey.GetVarCoeffs());
1156 
1157  MatrixKey invMasskey(
1160 
1161  DNekScalMat &invMass = *GetLocMatrix(invMasskey);
1162 
1163  Mat = Mat + Dmat * invMass * Transpose(Dmat);
1164  }
1165  }
1166  else if (mkey.HasVarCoeff(Coeffs[i]))
1167  {
1168  MatrixKey DmatkeyL(DerivType[i], DetShapeType(), *this,
1170  mkey.GetVarCoeffAsMap(Coeffs[i]));
1171 
1172  MatrixKey DmatkeyR(DerivType[i], DetShapeType(), *this);
1173 
1174  DNekScalMat &DmatL = *GetLocMatrix(DmatkeyL);
1175  DNekScalMat &DmatR = *GetLocMatrix(DmatkeyR);
1176  Mat = Mat + DmatL * invMass * Transpose(DmatR);
1177  }
1178  else
1179  {
1180  DNekScalMat &Dmat = *GetLocMatrix(DerivType[i]);
1181  Mat = Mat + Dmat * invMass * Transpose(Dmat);
1182  }
1183  }
1184 
1185  // Add Mass Matrix Contribution for Helmholtz problem
1187  Mat = Mat + lambdaval * Mass;
1188 
1189  // Add tau*E_l using elemental mass matrices on each edge
1190  for (i = 0; i < nedges; ++i)
1191  {
1192  EdgeExp = GetTraceExp(i);
1193  order_e = EdgeExp->GetNcoeffs();
1194 
1195  int nq = EdgeExp->GetNumPoints(0);
1196  GetTraceToElementMap(i, emap, sign, edgedir);
1197 
1198  // @TODO: Document
1199  StdRegions::VarCoeffMap edgeVarCoeffs;
1200  if (mkey.HasVarCoeff(StdRegions::eVarCoeffD00))
1201  {
1202  Array<OneD, NekDouble> mu(nq);
1204  i, EdgeExp, mkey.GetVarCoeff(StdRegions::eVarCoeffD00),
1205  mu);
1206  edgeVarCoeffs[StdRegions::eVarCoeffMass] = mu;
1207  }
1208  DNekScalMat &eMass = *EdgeExp->GetLocMatrix(
1210  edgeVarCoeffs);
1211  // DNekScalMat &eMass =
1212  // *EdgeExp->GetLocMatrix(StdRegions::eMass);
1213 
1214  for (j = 0; j < order_e; ++j)
1215  {
1216  for (k = 0; k < order_e; ++k)
1217  {
1218  Mat(emap[j], emap[k]) =
1219  Mat(emap[j], emap[k]) +
1220  tau * sign[j] * sign[k] * eMass(j, k);
1221  }
1222  }
1223  }
1224  }
1225  break;
1226  // U^e (P22)
1228  {
1229  int i, j, k;
1230  int nbndry = NumDGBndryCoeffs();
1231  int ncoeffs = GetNcoeffs();
1232  int nedges = GetNtraces();
1233  NekDouble tau = mkey.GetConstFactor(StdRegions::eFactorTau);
1234 
1235  Array<OneD, NekDouble> lambda(nbndry);
1236  DNekVec Lambda(nbndry, lambda, eWrapper);
1237  Array<OneD, NekDouble> ulam(ncoeffs);
1238  DNekVec Ulam(ncoeffs, ulam, eWrapper);
1239  Array<OneD, NekDouble> f(ncoeffs);
1240  DNekVec F(ncoeffs, f, eWrapper);
1241 
1242  Array<OneD, ExpansionSharedPtr> EdgeExp(nedges);
1243  // declare matrix space
1244  returnval =
1246  DNekMat &Umat = *returnval;
1247 
1248  // Z^e matrix
1249  MatrixKey newkey(StdRegions::eInvHybridDGHelmholtz, DetShapeType(),
1250  *this, mkey.GetConstFactors(),
1251  mkey.GetVarCoeffs());
1252  DNekScalMat &invHmat = *GetLocMatrix(newkey);
1253 
1254  Array<OneD, unsigned int> emap;
1255  Array<OneD, int> sign;
1256 
1257  for (i = 0; i < nedges; ++i)
1258  {
1259  EdgeExp[i] = GetTraceExp(i);
1260  }
1261 
1262  // for each degree of freedom of the lambda space
1263  // calculate Umat entry
1264  // Generate Lambda to U_lambda matrix
1265  for (j = 0; j < nbndry; ++j)
1266  {
1267  // standard basis vectors e_j
1268  Vmath::Zero(nbndry, &lambda[0], 1);
1269  Vmath::Zero(ncoeffs, &f[0], 1);
1270  lambda[j] = 1.0;
1271 
1272  SetTraceToGeomOrientation(EdgeExp, lambda);
1273 
1274  // Compute F = [I D_1 M^{-1} D_2 M^{-1}] C e_j
1275  AddHDGHelmholtzTraceTerms(tau, lambda, EdgeExp,
1276  mkey.GetVarCoeffs(), f);
1277 
1278  // Compute U^e_j
1279  Ulam = invHmat * F; // generate Ulam from lambda
1280 
1281  // fill column of matrix
1282  for (k = 0; k < ncoeffs; ++k)
1283  {
1284  Umat(k, j) = Ulam[k];
1285  }
1286  }
1287  }
1288  break;
1289  // Q_0, Q_1, Q_2 matrices (P23)
1290  // Each are a product of a row of Eqn 32 with the C matrix.
1291  // Rather than explicitly computing all of Eqn 32, we note each
1292  // row is almost a multiple of U^e, so use that as our starting
1293  // point.
1297  {
1298  int i = 0;
1299  int j = 0;
1300  int k = 0;
1301  int dir = 0;
1302  int nbndry = NumDGBndryCoeffs();
1303  int ncoeffs = GetNcoeffs();
1304  int nedges = GetNtraces();
1305  int shapedim = 2;
1306 
1307  Array<OneD, NekDouble> lambda(nbndry);
1308  DNekVec Lambda(nbndry, lambda, eWrapper);
1309  Array<OneD, ExpansionSharedPtr> EdgeExp(nedges);
1310 
1311  Array<OneD, NekDouble> ulam(ncoeffs);
1312  DNekVec Ulam(ncoeffs, ulam, eWrapper);
1313  Array<OneD, NekDouble> f(ncoeffs);
1314  DNekVec F(ncoeffs, f, eWrapper);
1315 
1316  // declare matrix space
1317  returnval =
1319  DNekMat &Qmat = *returnval;
1320 
1321  // Lambda to U matrix
1322  MatrixKey lamToUkey(StdRegions::eHybridDGLamToU, DetShapeType(),
1323  *this, mkey.GetConstFactors(),
1324  mkey.GetVarCoeffs());
1325  DNekScalMat &lamToU = *GetLocMatrix(lamToUkey);
1326 
1327  // Inverse mass matrix
1329 
1330  for (i = 0; i < nedges; ++i)
1331  {
1332  EdgeExp[i] = GetTraceExp(i);
1333  }
1334 
1335  // Weak Derivative matrix
1336  DNekScalMatSharedPtr Dmat;
1337  switch (mkey.GetMatrixType())
1338  {
1340  dir = 0;
1342  break;
1344  dir = 1;
1346  break;
1348  dir = 2;
1350  break;
1351  default:
1352  ASSERTL0(false, "Direction not known");
1353  break;
1354  }
1355 
1356  const StdRegions::VarCoeffMap &varcoeffs = mkey.GetVarCoeffs();
1357  if (varcoeffs.find(StdRegions::eVarCoeffMF1x) != varcoeffs.end())
1358  {
1359  StdRegions::VarCoeffMap VarCoeffDirDeriv;
1360  VarCoeffDirDeriv[StdRegions::eVarCoeffMF] =
1361  GetMF(dir, shapedim, varcoeffs);
1362  VarCoeffDirDeriv[StdRegions::eVarCoeffMFDiv] =
1363  GetMFDiv(dir, varcoeffs);
1364 
1365  MatrixKey Dmatkey(
1367  StdRegions::NullConstFactorMap, VarCoeffDirDeriv);
1368 
1369  Dmat = GetLocMatrix(Dmatkey);
1370 
1371  StdRegions::VarCoeffMap Weight;
1372  Weight[StdRegions::eVarCoeffMass] =
1373  GetMFMag(dir, mkey.GetVarCoeffs());
1374 
1375  MatrixKey invMasskey(StdRegions::eInvMass, DetShapeType(),
1377  Weight);
1378 
1379  invMass = *GetLocMatrix(invMasskey);
1380  }
1381  else
1382  {
1386 
1387  Dmat = GetLocMatrix(DerivType[dir]);
1388 
1389  MatrixKey invMasskey(StdRegions::eInvMass, DetShapeType(),
1390  *this);
1391  invMass = *GetLocMatrix(invMasskey);
1392  }
1393 
1394  // for each degree of freedom of the lambda space
1395  // calculate Qmat entry
1396  // Generate Lambda to Q_lambda matrix
1397  for (j = 0; j < nbndry; ++j)
1398  {
1399  Vmath::Zero(nbndry, &lambda[0], 1);
1400  lambda[j] = 1.0;
1401 
1402  // for lambda[j] = 1 this is the solution to ulam
1403  for (k = 0; k < ncoeffs; ++k)
1404  {
1405  Ulam[k] = lamToU(k, j);
1406  }
1407 
1408  // -D^T ulam
1409  Vmath::Neg(ncoeffs, &ulam[0], 1);
1410  F = Transpose(*Dmat) * Ulam;
1411 
1412  SetTraceToGeomOrientation(EdgeExp, lambda);
1413 
1414  // Add the C terms resulting from the I's on the
1415  // diagonals of Eqn 32
1416  AddNormTraceInt(dir, lambda, EdgeExp, f, mkey.GetVarCoeffs());
1417 
1418  // finally multiply by inverse mass matrix
1419  Ulam = invMass * F;
1420 
1421  // fill column of matrix (Qmat is in column major format)
1422  Vmath::Vcopy(ncoeffs, &ulam[0], 1,
1423  &(Qmat.GetPtr())[0] + j * ncoeffs, 1);
1424  }
1425  }
1426  break;
1427  // Matrix K (P23)
1429  {
1430  int i, j, e, cnt;
1431  int order_e, nquad_e;
1432  int nbndry = NumDGBndryCoeffs();
1433  int coordim = GetCoordim();
1434  int nedges = GetNtraces();
1435  NekDouble tau = mkey.GetConstFactor(StdRegions::eFactorTau);
1436  StdRegions::VarCoeffMap::const_iterator x;
1437  const StdRegions::VarCoeffMap &varcoeffs = mkey.GetVarCoeffs();
1438  bool mmf =
1439  (varcoeffs.find(StdRegions::eVarCoeffMF1x) != varcoeffs.end());
1440 
1441  Array<OneD, NekDouble> work, varcoeff_work;
1442  Array<OneD, const Array<OneD, NekDouble>> normals;
1443  Array<OneD, ExpansionSharedPtr> EdgeExp(nedges);
1444  Array<OneD, NekDouble> lam(nbndry);
1445 
1446  Array<OneD, unsigned int> emap;
1447  Array<OneD, int> sign;
1448  StdRegions::Orientation edgedir;
1449 
1450  // declare matrix space
1451  returnval =
1453  DNekMat &BndMat = *returnval;
1454 
1455  DNekScalMatSharedPtr LamToQ[3];
1456 
1457  // Matrix to map Lambda to U
1458  MatrixKey LamToUkey(StdRegions::eHybridDGLamToU, DetShapeType(),
1459  *this, mkey.GetConstFactors(),
1460  mkey.GetVarCoeffs());
1461  DNekScalMat &LamToU = *GetLocMatrix(LamToUkey);
1462 
1463  // Matrix to map Lambda to Q0
1464  MatrixKey LamToQ0key(StdRegions::eHybridDGLamToQ0, DetShapeType(),
1465  *this, mkey.GetConstFactors(),
1466  mkey.GetVarCoeffs());
1467  LamToQ[0] = GetLocMatrix(LamToQ0key);
1468 
1469  // Matrix to map Lambda to Q1
1470  MatrixKey LamToQ1key(StdRegions::eHybridDGLamToQ1, DetShapeType(),
1471  *this, mkey.GetConstFactors(),
1472  mkey.GetVarCoeffs());
1473  LamToQ[1] = GetLocMatrix(LamToQ1key);
1474 
1475  // Matrix to map Lambda to Q2 for 3D coordinates
1476  if (coordim == 3)
1477  {
1478  MatrixKey LamToQ2key(
1480  mkey.GetConstFactors(), mkey.GetVarCoeffs());
1481  LamToQ[2] = GetLocMatrix(LamToQ2key);
1482  }
1483 
1484  // Set up edge segment expansions from local geom info
1485  for (i = 0; i < nedges; ++i)
1486  {
1487  EdgeExp[i] = GetTraceExp(i);
1488  }
1489 
1490  // Set up matrix derived from <mu, Q_lam.n - \tau (U_lam - Lam) >
1491  for (i = 0; i < nbndry; ++i)
1492  {
1493  cnt = 0;
1494 
1495  Vmath::Zero(nbndry, lam, 1);
1496  lam[i] = 1.0;
1497  SetTraceToGeomOrientation(EdgeExp, lam);
1498 
1499  for (e = 0; e < nedges; ++e)
1500  {
1501  order_e = EdgeExp[e]->GetNcoeffs();
1502  nquad_e = EdgeExp[e]->GetNumPoints(0);
1503 
1504  normals = GetTraceNormal(e);
1505  edgedir = GetTraceOrient(e);
1506 
1507  work = Array<OneD, NekDouble>(nquad_e);
1508  varcoeff_work = Array<OneD, NekDouble>(nquad_e);
1509 
1510  GetTraceToElementMap(e, emap, sign, edgedir);
1511 
1512  StdRegions::VarCoeffType VarCoeff[3] = {
1515 
1516  // Q0 * n0 (BQ_0 terms)
1517  Array<OneD, NekDouble> edgeCoeffs(order_e);
1518  Array<OneD, NekDouble> edgePhys(nquad_e);
1519  for (j = 0; j < order_e; ++j)
1520  {
1521  edgeCoeffs[j] = sign[j] * (*LamToQ[0])(emap[j], i);
1522  }
1523 
1524  EdgeExp[e]->BwdTrans(edgeCoeffs, edgePhys);
1525  // @TODO Var coeffs
1526  // Multiply by variable coefficient
1527  // if ((x =
1528  // varcoeffs.find(VarCoeff[0]))
1529  // != varcoeffs.end())
1530  // {
1531  // GetPhysEdgeVarCoeffsFromElement(e,EdgeExp[e],x->second,varcoeff_work);
1532  // Vmath::Vmul(nquad_e,varcoeff_work,1,EdgeExp[e]->GetPhys(),1,EdgeExp[e]->UpdatePhys(),1);
1533  // }
1534  if (mmf)
1535  {
1536  Array<OneD, NekDouble> ncdotMF = GetnEdgecdotMF(
1537  0, e, EdgeExp[e], normals, varcoeffs);
1538  Vmath::Vmul(nquad_e, ncdotMF, 1, edgePhys, 1, work, 1);
1539  }
1540  else
1541  {
1542  Vmath::Vmul(nquad_e, normals[0], 1, edgePhys, 1, work,
1543  1);
1544  }
1545 
1546  // Q1 * n1 (BQ_1 terms)
1547  for (j = 0; j < order_e; ++j)
1548  {
1549  edgeCoeffs[j] = sign[j] * (*LamToQ[1])(emap[j], i);
1550  }
1551 
1552  EdgeExp[e]->BwdTrans(edgeCoeffs, edgePhys);
1553 
1554  // @TODO var coeffs
1555  // Multiply by variable coefficients
1556  // if ((x =
1557  // varcoeffs.find(VarCoeff[1]))
1558  // != varcoeffs.end())
1559  // {
1560  // GetPhysEdgeVarCoeffsFromElement(e,EdgeExp[e],x->second,varcoeff_work);
1561  // Vmath::Vmul(nquad_e,varcoeff_work,1,EdgeExp[e]->GetPhys(),1,EdgeExp[e]->UpdatePhys(),1);
1562  // }
1563 
1564  if (mmf)
1565  {
1566  Array<OneD, NekDouble> ncdotMF = GetnEdgecdotMF(
1567  1, e, EdgeExp[e], normals, varcoeffs);
1568  Vmath::Vvtvp(nquad_e, ncdotMF, 1, edgePhys, 1, work, 1,
1569  work, 1);
1570  }
1571  else
1572  {
1573  Vmath::Vvtvp(nquad_e, normals[1], 1, edgePhys, 1, work,
1574  1, work, 1);
1575  }
1576 
1577  // Q2 * n2 (BQ_2 terms)
1578  if (coordim == 3)
1579  {
1580  for (j = 0; j < order_e; ++j)
1581  {
1582  edgeCoeffs[j] = sign[j] * (*LamToQ[2])(emap[j], i);
1583  }
1584 
1585  EdgeExp[e]->BwdTrans(edgeCoeffs, edgePhys);
1586  // @TODO var coeffs
1587  // Multiply by variable coefficients
1588  // if ((x =
1589  // varcoeffs.find(VarCoeff[2]))
1590  // != varcoeffs.end())
1591  // {
1592  // GetPhysEdgeVarCoeffsFromElement(e,EdgeExp[e],x->second,varcoeff_work);
1593  // Vmath::Vmul(nquad_e,varcoeff_work,1,EdgeExp[e]->GetPhys(),1,EdgeExp[e]->UpdatePhys(),1);
1594  // }
1595 
1596  Vmath::Vvtvp(nquad_e, normals[2], 1, edgePhys, 1, work,
1597  1, work, 1);
1598  }
1599 
1600  // - tau (ulam - lam)
1601  // Corresponds to the G and BU terms.
1602  for (j = 0; j < order_e; ++j)
1603  {
1604  edgeCoeffs[j] =
1605  sign[j] * LamToU(emap[j], i) - lam[cnt + j];
1606  }
1607 
1608  EdgeExp[e]->BwdTrans(edgeCoeffs, edgePhys);
1609 
1610  // Multiply by variable coefficients
1611  if ((x = varcoeffs.find(VarCoeff[0])) != varcoeffs.end())
1612  {
1614  e, EdgeExp[e], x->second.GetValue(), varcoeff_work);
1615  Vmath::Vmul(nquad_e, varcoeff_work, 1, edgePhys, 1,
1616  edgePhys, 1);
1617  }
1618 
1619  Vmath::Svtvp(nquad_e, -tau, edgePhys, 1, work, 1, work, 1);
1620  /// TODO: Add variable coeffs
1621  EdgeExp[e]->IProductWRTBase(work, edgeCoeffs);
1622 
1623  EdgeExp[e]->SetCoeffsToOrientation(edgedir, edgeCoeffs,
1624  edgeCoeffs);
1625 
1626  for (j = 0; j < order_e; ++j)
1627  {
1628  BndMat(cnt + j, i) = edgeCoeffs[j];
1629  }
1630 
1631  cnt += order_e;
1632  }
1633  }
1634  }
1635  break;
1636  // HDG postprocessing
1638  {
1639  MatrixKey lapkey(StdRegions::eLaplacian, DetShapeType(), *this,
1640  mkey.GetConstFactors(), mkey.GetVarCoeffs());
1641  DNekScalMat &LapMat = *GetLocMatrix(lapkey);
1642 
1644  LapMat.GetRows(), LapMat.GetColumns());
1645  DNekMatSharedPtr lmat = returnval;
1646 
1647  (*lmat) = LapMat;
1648 
1649  // replace first column with inner product wrt 1
1650  int nq = GetTotPoints();
1651  Array<OneD, NekDouble> tmp(nq);
1652  Array<OneD, NekDouble> outarray(m_ncoeffs);
1653  Vmath::Fill(nq, 1.0, tmp, 1);
1654  IProductWRTBase(tmp, outarray);
1655 
1656  Vmath::Vcopy(m_ncoeffs, &outarray[0], 1, &(lmat->GetPtr())[0], 1);
1657  lmat->Invert();
1658  }
1659  break;
1661  {
1662  int ntraces = GetNtraces();
1663  int ncoords = GetCoordim();
1664  int nphys = GetTotPoints();
1665  Array<OneD, const Array<OneD, NekDouble>> normals;
1666  Array<OneD, NekDouble> phys(nphys);
1667  returnval =
1669  DNekMat &Mat = *returnval;
1670  Vmath::Zero(m_ncoeffs * m_ncoeffs, Mat.GetPtr(), 1);
1671 
1672  Array<OneD, Array<OneD, NekDouble>> Deriv(3, NullNekDouble1DArray);
1673 
1674  for (int d = 0; d < ncoords; ++d)
1675  {
1676  Deriv[d] = Array<OneD, NekDouble>(nphys);
1677  }
1678 
1679  Array<OneD, int> tracepts(ntraces);
1680  Array<OneD, ExpansionSharedPtr> traceExp(ntraces);
1681  int maxtpts = 0;
1682  for (int t = 0; t < ntraces; ++t)
1683  {
1684  traceExp[t] = GetTraceExp(t);
1685  tracepts[t] = traceExp[t]->GetTotPoints();
1686  maxtpts = (maxtpts > tracepts[t]) ? maxtpts : tracepts[t];
1687  }
1688 
1689  Array<OneD, NekDouble> val(maxtpts), tmp, tmp1;
1690 
1691  Array<OneD, Array<OneD, NekDouble>> dphidn(ntraces);
1692  for (int t = 0; t < ntraces; ++t)
1693  {
1694  dphidn[t] =
1695  Array<OneD, NekDouble>(m_ncoeffs * tracepts[t], 0.0);
1696  }
1697 
1698  for (int i = 0; i < m_ncoeffs; ++i)
1699  {
1700  FillMode(i, phys);
1701  PhysDeriv(phys, Deriv[0], Deriv[1], Deriv[2]);
1702 
1703  for (int t = 0; t < ntraces; ++t)
1704  {
1705  const NormalVector norm = GetTraceNormal(t);
1706 
1707  LibUtilities::BasisKey fromkey = GetTraceBasisKey(t);
1708  LibUtilities::BasisKey tokey =
1709  traceExp[t]->GetBasis(0)->GetBasisKey();
1710  bool DoInterp = (fromkey != tokey);
1711 
1712  Array<OneD, NekDouble> n(tracepts[t]);
1713  ;
1714  for (int d = 0; d < ncoords; ++d)
1715  {
1716  // if variable p may need to interpolate
1717  if (DoInterp)
1718  {
1719  LibUtilities::Interp1D(fromkey, norm[d], tokey, n);
1720  }
1721  else
1722  {
1723  n = norm[d];
1724  }
1725 
1726  GetTracePhysVals(t, traceExp[t], Deriv[d], val,
1727  v_GetTraceOrient(t));
1728 
1729  Vmath::Vvtvp(tracepts[t], n, 1, val, 1,
1730  tmp = dphidn[t] + i * tracepts[t], 1,
1731  tmp1 = dphidn[t] + i * tracepts[t], 1);
1732  }
1733  }
1734  }
1735 
1736  for (int t = 0; t < ntraces; ++t)
1737  {
1738  int nt = tracepts[t];
1739  NekDouble h, p;
1740  TraceNormLen(t, h, p);
1741 
1742  // scaling from GJP paper
1743  NekDouble scale =
1744  (p == 1) ? 0.02 * h * h : 0.8 * pow(p + 1, -4.0) * h * h;
1745 
1746  for (int i = 0; i < m_ncoeffs; ++i)
1747  {
1748  for (int j = i; j < m_ncoeffs; ++j)
1749  {
1750  Vmath::Vmul(nt, dphidn[t] + i * nt, 1,
1751  dphidn[t] + j * nt, 1, val, 1);
1752  Mat(i, j) =
1753  Mat(i, j) + scale * traceExp[t]->Integral(val);
1754  }
1755  }
1756  }
1757 
1758  // fill in symmetric components.
1759  for (int i = 0; i < m_ncoeffs; ++i)
1760  {
1761  for (int j = 0; j < i; ++j)
1762  {
1763  Mat(i, j) = Mat(j, i);
1764  }
1765  }
1766  }
1767  break;
1768  default:
1769  ASSERTL0(false,
1770  "This matrix type cannot be generated from this class");
1771  break;
1772  }
1773 
1774  return returnval;
1775 }
void SetTraceToGeomOrientation(Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, NekDouble > &inout)
void AddHDGHelmholtzTraceTerms(const NekDouble tau, const Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, const StdRegions::VarCoeffMap &dirForcing, Array< OneD, NekDouble > &outarray)
void GetTracePhysVals(const int trace, const StdRegions::StdExpansionSharedPtr &TraceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=StdRegions::eNoOrientation)
Definition: Expansion.h:202
ExpansionSharedPtr GetTraceExp(const int traceid)
Definition: Expansion.h:416
void TraceNormLen(const int traceid, NekDouble &h, NekDouble &p)
Definition: Expansion.h:252
void FillMode(const int mode, Array< OneD, NekDouble > &outarray)
This function fills the array outarray with the mode-th mode of the expansion.
Definition: StdExpansion.h:497
const LibUtilities::BasisKey GetTraceBasisKey(const int i, int k=-1) const
This function returns the basis key belonging to the i-th trace.
Definition: StdExpansion.h:305
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 expa...
Definition: StdExpansion.h:534
void PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Definition: StdExpansion.h:848
void Interp1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
this function interpolates a 1D function evaluated at the quadrature points of the basis fbasis0 to ...
Definition: Interp.cpp:52
Array< OneD, Array< OneD, NekDouble > > NormalVector
Definition: Expansion.h:53
static Array< OneD, NekDouble > NullNekDouble1DArray
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45

References AddHDGHelmholtzTraceTerms(), AddNormTraceInt(), Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, ASSERTL1, Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eFactorTau, Nektar::StdRegions::eForwards, Nektar::StdRegions::eHybridDGHelmBndLam, Nektar::StdRegions::eHybridDGHelmholtz, Nektar::StdRegions::eHybridDGLamToQ0, Nektar::StdRegions::eHybridDGLamToQ1, Nektar::StdRegions::eHybridDGLamToQ2, Nektar::StdRegions::eHybridDGLamToU, Nektar::StdRegions::eInvHybridDGHelmholtz, Nektar::StdRegions::eInvLaplacianWithUnityMean, Nektar::StdRegions::eInvMass, Nektar::StdRegions::eLaplacian, Nektar::StdRegions::eMass, Nektar::StdRegions::eNormDerivOnTrace, Nektar::StdRegions::eVarCoeffD00, Nektar::StdRegions::eVarCoeffD11, Nektar::StdRegions::eVarCoeffD22, Nektar::StdRegions::eVarCoeffMass, Nektar::StdRegions::eVarCoeffMF, Nektar::StdRegions::eVarCoeffMF1x, Nektar::StdRegions::eVarCoeffMFDiv, Nektar::StdRegions::eWeakDeriv0, Nektar::StdRegions::eWeakDeriv1, Nektar::StdRegions::eWeakDeriv2, Nektar::StdRegions::eWeakDirectionalDeriv, Nektar::eWrapper, Vmath::Fill(), Nektar::StdRegions::StdExpansion::FillMode(), Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdMatrixKey::GetConstFactors(), Nektar::StdRegions::StdExpansion::GetCoordim(), Nektar::LocalRegions::Expansion::GetLocMatrix(), Nektar::StdRegions::StdMatrixKey::GetMatrixType(), Nektar::LocalRegions::Expansion::GetMF(), Nektar::LocalRegions::Expansion::GetMFDiv(), Nektar::LocalRegions::Expansion::GetMFMag(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), GetnEdgecdotMF(), Nektar::StdRegions::StdExpansion::GetNtraces(), GetPhysEdgeVarCoeffsFromElement(), Nektar::StdRegions::StdExpansion::GetTotPoints(), Nektar::StdRegions::StdExpansion::GetTraceBasisKey(), Nektar::LocalRegions::Expansion::GetTraceExp(), Nektar::LocalRegions::Expansion::GetTraceNormal(), Nektar::LocalRegions::Expansion::GetTraceOrient(), Nektar::LocalRegions::Expansion::GetTracePhysVals(), Nektar::StdRegions::StdExpansion::GetTraceToElementMap(), Nektar::StdRegions::StdMatrixKey::GetVarCoeff(), Nektar::StdRegions::StdMatrixKey::GetVarCoeffAsMap(), Nektar::StdRegions::StdMatrixKey::GetVarCoeffs(), Nektar::StdRegions::StdMatrixKey::HasVarCoeff(), Nektar::LibUtilities::Interp1D(), Nektar::StdRegions::StdExpansion::IProductWRTBase(), Nektar::StdRegions::StdExpansion::IsBoundaryInteriorExpansion(), Nektar::StdRegions::StdExpansion::m_ncoeffs, Vmath::Neg(), Nektar::StdRegions::NullConstFactorMap, Nektar::NullNekDouble1DArray, Nektar::StdRegions::StdExpansion::NumDGBndryCoeffs(), CellMLToNektar.cellml_metadata::p, Nektar::StdRegions::StdExpansion::PhysDeriv(), SetTraceToGeomOrientation(), sign, Vmath::Svtvp(), Nektar::LocalRegions::Expansion::TraceNormLen(), Nektar::Transpose(), Nektar::LocalRegions::Expansion::v_GetTraceOrient(), Vmath::Vcopy(), Vmath::Vmul(), Vmath::Vvtvp(), and Vmath::Zero().

Referenced by CreateMatrix(), Nektar::LocalRegions::NodalTriExp::v_GenMatrix(), Nektar::LocalRegions::QuadExp::v_GenMatrix(), and Nektar::LocalRegions::TriExp::v_GenMatrix().

◆ v_GenTraceExp()

void Nektar::LocalRegions::Expansion2D::v_GenTraceExp ( const int  traceid,
ExpansionSharedPtr exp 
)
overridevirtual

◆ v_ReOrientTracePhysMap()

void Nektar::LocalRegions::Expansion2D::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 2064 of file Expansion2D.cpp.

2067 {
2068  boost::ignore_unused(nq1);
2069 
2070  if (idmap.size() != nq0)
2071  {
2072  idmap = Array<OneD, int>(nq0);
2073  }
2074  switch (orient)
2075  {
2076  case StdRegions::eForwards:
2077  // Fwd
2078  for (int i = 0; i < nq0; ++i)
2079  {
2080  idmap[i] = i;
2081  }
2082  break;
2084  {
2085  // Bwd
2086  for (int i = 0; i < nq0; ++i)
2087  {
2088  idmap[i] = nq0 - 1 - i;
2089  }
2090  }
2091  break;
2092  default:
2093  ASSERTL0(false, "Unknown orientation");
2094  break;
2095  }
2096 }

References ASSERTL0, Nektar::StdRegions::eBackwards, and Nektar::StdRegions::eForwards.

◆ v_SetUpPhysNormals()

void Nektar::LocalRegions::Expansion2D::v_SetUpPhysNormals ( const int  edge)
overrideprotectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 2059 of file Expansion2D.cpp.

2060 {
2061  v_ComputeTraceNormal(edge);
2062 }
virtual void v_ComputeTraceNormal(const int id)
Definition: Expansion.cpp:872

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

◆ v_TraceNormLen()

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

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 2154 of file Expansion2D.cpp.

2155 {
2157 
2158  int nverts = geom->GetNumVerts();
2159 
2160  // vertices on edges
2161  SpatialDomains::PointGeom ev0 = *geom->GetVertex(traceid);
2162  SpatialDomains::PointGeom ev1 = *geom->GetVertex((traceid + 1) % nverts);
2163 
2164  // vertex on adjacent edge to ev0
2165  SpatialDomains::PointGeom vadj =
2166  *geom->GetVertex((traceid + (nverts - 1)) % nverts);
2167 
2168  // calculate perpendicular distance of normal length
2169  // from first vertex
2170  NekDouble h1 = ev0.dist(vadj);
2171  SpatialDomains::PointGeom Dx, Dx1;
2172 
2173  Dx.Sub(ev1, ev0);
2174  Dx1.Sub(vadj, ev0);
2175 
2176  NekDouble d1 = Dx.dot(Dx1);
2177  NekDouble lenDx = Dx.dot(Dx);
2178  h = sqrt(h1 * h1 - d1 * d1 / lenDx);
2179 
2180  // perpendicular distanace from second vertex
2181  SpatialDomains::PointGeom vadj1 = *geom->GetVertex((traceid + 2) % nverts);
2182 
2183  h1 = ev1.dist(vadj1);
2184  Dx1.Sub(vadj1, ev1);
2185  d1 = Dx.dot(Dx1);
2186 
2187  h = (h + sqrt(h1 * h1 - d1 * d1 / lenDx)) * 0.5;
2188 
2189  int dirn = (geom->GetDir(traceid) == 0) ? 1 : 0;
2190 
2191  p = (NekDouble)(GetBasisNumModes(dirn) - 1);
2192 }
std::shared_ptr< Geometry > GeometrySharedPtr
Definition: Geometry.h:53
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References Nektar::SpatialDomains::PointGeom::dist(), Nektar::SpatialDomains::PointGeom::dot(), Nektar::StdRegions::StdExpansion::GetBasisNumModes(), Nektar::LocalRegions::Expansion::GetGeom(), Nektar::SpatialDomains::Geometry::GetVertex(), CellMLToNektar.cellml_metadata::p, tinysimd::sqrt(), and Nektar::SpatialDomains::PointGeom::Sub().

◆ v_VectorFlux()

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

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 2138 of file Expansion2D.cpp.

2140 {
2141  const Array<OneD, const Array<OneD, NekDouble>> &normals =
2142  GetLeftAdjacentElementExp()->GetTraceNormal(
2144 
2145  int nq = GetTotPoints();
2146  Array<OneD, NekDouble> Fn(nq);
2147  Vmath::Vmul(nq, &vec[0][0], 1, &normals[0][0], 1, &Fn[0], 1);
2148  Vmath::Vvtvp(nq, &vec[1][0], 1, &normals[1][0], 1, &Fn[0], 1, &Fn[0], 1);
2149  Vmath::Vvtvp(nq, &vec[2][0], 1, &normals[2][0], 1, &Fn[0], 1, &Fn[0], 1);
2150 
2151  return StdExpansion::Integral(Fn);
2152 }
ExpansionSharedPtr GetLeftAdjacentElementExp() const
Definition: Expansion.h:443
int GetLeftAdjacentElementTrace() const
Definition: Expansion.h:456

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

Member Data Documentation

◆ m_requireNeg

std::vector<bool> Nektar::LocalRegions::Expansion2D::m_requireNeg
protected

Definition at line 118 of file Expansion2D.h.

Referenced by v_AddEdgeNormBoundaryInt().