Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
Nektar::LocalRegions::Expansion2D Class Reference

#include <Expansion2D.h>

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

Public Member Functions

 Expansion2D (SpatialDomains::Geometry2DSharedPtr pGeom)
virtual ~Expansion2D ()
void SetTraceToGeomOrientation (Array< OneD, StdRegions::StdExpansionSharedPtr > &EdgeExp, Array< OneD, NekDouble > &inout)
ExpansionSharedPtr GetEdgeExp (int edge, bool SetUpNormal=true)
void SetEdgeExp (const int edge, ExpansionSharedPtr &e)
void AddNormTraceInt (const int dir, Array< OneD, StdRegions::StdExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &edgeCoeffs, Array< OneD, NekDouble > &outarray)
void AddNormTraceInt (const int dir, Array< OneD, const NekDouble > &inarray, Array< OneD, StdRegions::StdExpansionSharedPtr > &EdgeExp, Array< OneD, NekDouble > &outarray, const StdRegions::VarCoeffMap &varcoeffs)
void AddEdgeBoundaryInt (const int edge, StdRegions::StdExpansionSharedPtr &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, StdRegions::StdExpansionSharedPtr > &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, StdRegions::StdExpansionSharedPtr > &EdgeExp, const StdRegions::VarCoeffMap &dirForcing, Array< OneD, NekDouble > &outarray)
Expansion3DSharedPtr GetLeftAdjacentElementExp () const
Expansion3DSharedPtr GetRightAdjacentElementExp () const
int GetLeftAdjacentElementFace () const
int GetRightAdjacentElementFace () const
void SetAdjacentElementExp (int face, Expansion3DSharedPtr &f)
SpatialDomains::Geometry2DSharedPtr GetGeom2D () const
- Public Member Functions inherited from Nektar::LocalRegions::Expansion
 Expansion (SpatialDomains::GeometrySharedPtr pGeom)
 Expansion (const Expansion &pSrc)
virtual ~Expansion ()
DNekScalMatSharedPtr GetLocMatrix (const LocalRegions::MatrixKey &mkey)
DNekScalMatSharedPtr GetLocMatrix (const StdRegions::MatrixType mtype, const StdRegions::ConstFactorMap &factors=StdRegions::NullConstFactorMap, const StdRegions::VarCoeffMap &varcoeffs=StdRegions::NullVarCoeffMap)
SpatialDomains::GeometrySharedPtr GetGeom () const
virtual const
SpatialDomains::GeomFactorsSharedPtr
v_GetMetricInfo () const
DNekMatSharedPtr BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
DNekMatSharedPtr BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd)
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion
 StdExpansion ()
 Default Constructor.
 StdExpansion (const int numcoeffs, const int numbases, const LibUtilities::BasisKey &Ba=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bb=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bc=LibUtilities::NullBasisKey)
 Constructor.
 StdExpansion (const StdExpansion &T)
 Copy Constructor.
virtual ~StdExpansion ()
 Destructor.
int GetNumBases () const
 This function returns the number of 1D bases used in the expansion.
const Array< OneD, const
LibUtilities::BasisSharedPtr > & 
GetBase () const
 This function gets the shared point to basis.
const
LibUtilities::BasisSharedPtr
GetBasis (int dir) const
 This function gets the shared point to basis in the dir direction.
int GetNcoeffs (void) const
 This function returns the total number of coefficients used in the expansion.
int GetTotPoints () const
 This function returns the total number of quadrature points used in the element.
LibUtilities::BasisType GetBasisType (const int dir) const
 This function returns the type of basis used in the dir direction.
int GetBasisNumModes (const int dir) const
 This function returns the number of expansion modes in the dir direction.
int EvalBasisNumModesMax (void) const
 This function returns the maximum number of expansion modes over all local directions.
LibUtilities::PointsType GetPointsType (const int dir) const
 This function returns the type of quadrature points used in the dir direction.
int GetNumPoints (const int dir) const
 This function returns the number of quadrature points in the dir direction.
const Array< OneD, const
NekDouble > & 
GetPoints (const int dir) const
 This function returns a pointer to the array containing the quadrature points in dir direction.
int GetNverts () const
 This function returns the number of vertices of the expansion domain.
int GetNedges () const
 This function returns the number of edges of the expansion domain.
int GetEdgeNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th edge.
int GetTotalEdgeIntNcoeffs () const
int GetEdgeNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th edge.
int DetCartesianDirOfEdge (const int edge)
const LibUtilities::BasisKey DetEdgeBasisKey (const int i) const
const LibUtilities::BasisKey DetFaceBasisKey (const int i, const int k) const
int GetFaceNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th face.
int GetFaceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th face.
int GetFaceIntNcoeffs (const int i) const
int GetTotalFaceIntNcoeffs () const
LibUtilities::PointsKey GetFacePointsKey (const int i, const int j) const
int NumBndryCoeffs (void) const
int NumDGBndryCoeffs (void) const
LibUtilities::BasisType GetEdgeBasisType (const int i) const
 This function returns the type of expansion basis on the i-th edge.
int GetNfaces () const
 This function returns the number of faces of the expansion domain.
LibUtilities::ShapeType DetShapeType () const
 This function returns the shape of the expansion domain.
int GetShapeDimension () const
bool IsBoundaryInteriorExpansion ()
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space.
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space.
void FwdTrans_BndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
NekDouble Integral (const Array< OneD, const NekDouble > &inarray)
 This function integrates the specified function over the domain.
void FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 This function fills the array outarray with the mode-th mode of the expansion.
void IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 this function calculates the inner product of a given function f with the different modes of the expansion
void IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check)
void IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
int GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id.
void SetElmtId (const int id)
 Set the element id of this expansion when used in a list by returning value of m_elmt_id.
void GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2=NullNekDouble1DArray, Array< OneD, NekDouble > &coords_3=NullNekDouble1DArray)
 this function returns the physical coordinates of the quadrature points of the expansion
void GetCoord (const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
 given the coordinates of a point of the element in the local collapsed coordinate system, this function calculates the physical coordinates of the point
DNekMatSharedPtr GetStdMatrix (const StdMatrixKey &mkey)
DNekBlkMatSharedPtr GetStdStaticCondMatrix (const StdMatrixKey &mkey)
IndexMapValuesSharedPtr GetIndexMap (const IndexMapKey &ikey)
const Array< OneD, const
NekDouble > & 
GetPhysNormals (void)
void SetPhysNormals (Array< OneD, const NekDouble > &normal)
virtual void SetUpPhysNormals (const int edge)
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const 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)
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
void DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
StdRegions::Orientation GetFaceOrient (int face)
StdRegions::Orientation GetEorient (int edge)
StdRegions::Orientation GetPorient (int point)
StdRegions::Orientation GetCartesianEorient (int edge)
void SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
void SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
int CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
void ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs)
NekDouble StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
void AddEdgeNormBoundaryInt (const int edge, boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
void AddEdgeNormBoundaryInt (const int edge, boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
void AddEdgeNormBoundaryBiInt (const int edge, boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
void AddFaceNormBoundaryInt (const int face, boost::shared_ptr< StdExpansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
int GetCoordim ()
void GetBoundaryMap (Array< OneD, unsigned int > &outarray)
void GetInteriorMap (Array< OneD, unsigned int > &outarray)
int GetVertexMap (const int localVertexId, bool useCoeffPacking=false)
void GetEdgeInteriorMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
void GetFaceInteriorMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
void GetEdgeToElementMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
void GetFaceToElementMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int nummodesA=-1, int nummodesB=-1)
void GetEdgePhysVals (const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Extract the physical values along edge edge from inarray into outarray following the local edge orientation and point distribution defined by defined in EdgeExp.
void GetEdgePhysVals (const int edge, const boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void GetTracePhysVals (const int edge, const boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void GetVertexPhysVals (const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray)
void GetEdgeInterpVals (const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void GetEdgeQFactors (const int edge, Array< OneD, NekDouble > &outarray)
 Extract the metric factors to compute the contravariant fluxes along edge edge and stores them into outarray following the local edge orientation (i.e. anticlockwise convention).
void GetFacePhysVals (const int face, const boost::shared_ptr< StdExpansion > &FaceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=eNoOrientation)
void MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
DNekMatSharedPtr CreateGeneralMatrix (const StdMatrixKey &mkey)
 this function generates the mass matrix $\mathbf{M}[i][j] = \int \phi_i(\mathbf{x}) \phi_j(\mathbf{x}) d\mathbf{x}$
void GeneralMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
void LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void WeakDirectionalDerivMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void MassLevelCurvatureMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LinearAdvectionDiffusionReactionMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
void HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
DNekMatSharedPtr GenMatrix (const StdMatrixKey &mkey)
void PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void PhysDeriv_s (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds)
void PhysDeriv_n (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn)
void PhysDirectionalDeriv (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &outarray)
void StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
void StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void AddRobinMassMatrix (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
void AddRobinEdgeContribution (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, Array< OneD, NekDouble > &coeffs)
void DGDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, boost::shared_ptr< StdExpansion > > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &coeffs, 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.
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.
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta.
const boost::shared_ptr
< SpatialDomains::GeomFactors > & 
GetMetricInfo (void) const
virtual int v_GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id.
virtual const Array< OneD,
const NekDouble > & 
v_GetPhysNormals (void)
virtual void v_SetPhysNormals (Array< OneD, const NekDouble > &normal)
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
virtual void v_ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmode_offset, NekDouble *coeffs)
 Unpack data from input file assuming it comes from the same expansion type.
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 DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
virtual void v_DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
virtual StdRegions::Orientation v_GetFaceOrient (int face)
virtual StdRegions::Orientation v_GetEorient (int edge)
virtual StdRegions::Orientation v_GetCartesianEorient (int edge)
virtual StdRegions::Orientation v_GetPorient (int point)
NekDouble Linf (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete $ L_\infty$ error $ |\epsilon|_\infty = \max |u - u_{exact}|$ where $ u_{exact}$ is given by the array sol.
NekDouble L2 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete $ L_2$ error, $ | \epsilon |_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 dx \right]^{1/2} d\xi_1 $ where $ u_{exact}$ is given by the array sol.
NekDouble H1 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete $ H^1$ error, $ | \epsilon |^1_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 + \nabla(u - u_{exact})\cdot\nabla(u - u_{exact})\cdot dx \right]^{1/2} d\xi_1 $ where $ u_{exact}$ is given by the array sol.
const NormalVectorGetEdgeNormal (const int edge) const
void ComputeEdgeNormal (const int edge)
void NegateEdgeNormal (const int edge)
bool EdgeNormalNegated (const int edge)
void ComputeFaceNormal (const int face)
void NegateFaceNormal (const int face)
void ComputeVertexNormal (const int vertex)
const NormalVectorGetFaceNormal (const int face) const
const NormalVectorGetVertexNormal (const int vertex) const
const NormalVectorGetSurfaceNormal (const int id) const
const LibUtilities::PointsKeyVector GetPointsKeys () const
Array< OneD, unsigned int > GetEdgeInverseBoundaryMap (int eid)
Array< OneD, unsigned int > GetFaceInverseBoundaryMap (int fid, StdRegions::Orientation faceOrient=eNoOrientation)
DNekMatSharedPtr BuildInverseTransformationMatrix (const DNekScalMatSharedPtr &m_transformationmatrix)
template<class T >
boost::shared_ptr< T > as ()
- 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 ()
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.
NekDouble Integral (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &w0, const Array< OneD, const NekDouble > &w1)
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 DNekMatSharedPtr v_GenMatrix (const StdRegions::StdMatrixKey &mkey)
virtual void v_DGDeriv (const int dir, const Array< OneD, const NekDouble > &incoeffs, Array< OneD, StdRegions::StdExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &edgeCoeffs, Array< OneD, NekDouble > &out_d)
virtual void v_AddEdgeNormBoundaryInt (const int edge, StdRegions::StdExpansionSharedPtr &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
virtual void v_AddEdgeNormBoundaryInt (const int edge, StdRegions::StdExpansionSharedPtr &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
virtual void v_AddRobinMassMatrix (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
virtual void v_AddRobinEdgeContribution (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, Array< OneD, NekDouble > &coeffs)
virtual DNekMatSharedPtr v_BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd)
void GetPhysEdgeVarCoeffsFromElement (const int edge, StdRegions::StdExpansionSharedPtr &EdgeExp, const Array< OneD, const NekDouble > &varcoeff, Array< OneD, NekDouble > &outarray)
Array< OneD, unsigned int > v_GetEdgeInverseBoundaryMap (int eid)
virtual void v_NegateEdgeNormal (const int edge)
virtual bool v_EdgeNormalNegated (const int edge)
virtual void v_SetUpPhysNormals (const int edge)
const StdRegions::NormalVectorv_GetEdgeNormal (const int edge) const
const StdRegions::NormalVectorv_GetSurfaceNormal (const int id) const
- Protected Member Functions inherited from Nektar::LocalRegions::Expansion
void ComputeLaplacianMetric ()
void ComputeQuadratureMetric ()
virtual void v_MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_ComputeLaplacianMetric ()
virtual void v_GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3)
virtual DNekScalMatSharedPtr v_GetLocMatrix (const LocalRegions::MatrixKey &mkey)
virtual DNekMatSharedPtr v_BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
- 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.
IndexMapValuesSharedPtr CreateIndexMap (const IndexMapKey &ikey)
 Create an IndexMap which contains mapping information linking any specific element shape with either its boundaries, edges, faces, verteces, etc.
void BwdTrans_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void IProductWRTBase_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 GeneralMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void MassMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
void LaplacianMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LaplacianMatrixOp_MatFree (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void WeakDerivMatrixOp_MatFree (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void WeakDirectionalDerivMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void MassLevelCurvatureMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LinearAdvectionDiffusionReactionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
void HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void HelmholtzMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
virtual NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
virtual void v_AddEdgeNormBoundaryInt (const int edge, boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
virtual void v_AddEdgeNormBoundaryInt (const int edge, boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
virtual void v_AddFaceNormBoundaryInt (const int face, boost::shared_ptr< StdExpansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion2D
virtual NekDouble v_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.
virtual NekDouble v_PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
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)
virtual void v_HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)

Protected Attributes

std::vector< ExpansionWeakPtrm_edgeExp
std::vector< bool > m_requireNeg
std::map< int,
StdRegions::NormalVector
m_edgeNormals
std::map< int, bool > m_negatedNormals
Expansion3DWeakPtr m_elementLeft
Expansion3DWeakPtr m_elementRight
int m_elementFaceLeft
int m_elementFaceRight
- Protected Attributes inherited from Nektar::LocalRegions::Expansion
SpatialDomains::GeometrySharedPtr m_geom
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
MetricMap m_metrics
- Protected Attributes inherited from Nektar::StdRegions::StdExpansion
Array< OneD,
LibUtilities::BasisSharedPtr
m_base
int m_elmt_id
int m_ncoeffs
LibUtilities::NekManager
< StdMatrixKey, DNekMat,
StdMatrixKey::opLess
m_stdMatrixManager
LibUtilities::NekManager
< StdMatrixKey, DNekBlkMat,
StdMatrixKey::opLess
m_stdStaticCondMatrixManager
LibUtilities::NekManager
< IndexMapKey, IndexMapValues,
IndexMapKey::opLess
m_IndexMapManager

Detailed Description

Definition at line 58 of file Expansion2D.h.

Constructor & Destructor Documentation

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

Definition at line 47 of file Expansion2D.cpp.

References m_elementFaceLeft, and m_elementFaceRight.

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

Definition at line 62 of file Expansion2D.h.

{}

Member Function Documentation

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

For a given edge add the {F}_1j contributions

Variable coeffs

Definition at line 417 of file Expansion2D.cpp.

References Nektar::StdRegions::eVarCoeffD00, Nektar::StdRegions::eVarCoeffD11, Nektar::StdRegions::eVarCoeffD22, Nektar::StdRegions::StdExpansion::GetEdgeToElementMap(), GetPhysEdgeVarCoeffsFromElement(), sign, Nektar::StdRegions::StdExpansion::v_GetEorient(), and Vmath::Vmul().

Referenced by AddNormTraceInt().

{
int i;
int order_e = EdgeExp->GetNcoeffs();
int nquad_e = EdgeExp->GetNumPoints(0);
Array<OneD,unsigned int> map;
Array<OneD,int> sign;
Array<OneD, NekDouble> coeff(order_e);
GetEdgeToElementMap(edge, v_GetEorient(edge), map, sign);
StdRegions::VarCoeffMap::const_iterator x;
/// @TODO Variable coeffs
if ((x = varcoeffs.find(VarCoeff[0])) != varcoeffs.end())
{
Array<OneD, NekDouble> work(nquad_e);
edge, EdgeExp, x->second, work);
Vmath::Vmul(nquad_e, work, 1, edgePhys, 1, edgePhys, 1);
}
EdgeExp->IProductWRTBase(edgePhys, coeff);
// add data to out array
for(i = 0; i < order_e; ++i)
{
outarray[map[i]] += sign[i]*coeff[i];
}
}
void Nektar::LocalRegions::Expansion2D::AddHDGHelmholtzEdgeTerms ( const NekDouble  tau,
const int  edge,
Array< OneD, StdRegions::StdExpansionSharedPtr > &  EdgeExp,
Array< OneD, NekDouble > &  edgePhys,
const StdRegions::VarCoeffMap dirForcing,
Array< OneD, NekDouble > &  outarray 
)
inline

: What direction to use here??

: Document this (probably not needed)

Definition at line 491 of file Expansion2D.cpp.

References Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::StdRegions::eInvMass, Nektar::StdRegions::eVarCoeffD00, Nektar::StdRegions::eVarCoeffD11, Nektar::StdRegions::eVarCoeffD22, Nektar::StdRegions::eWeakDeriv0, Nektar::StdRegions::eWeakDeriv1, Nektar::StdRegions::eWeakDeriv2, Nektar::eWrapper, Nektar::StdRegions::StdExpansion::GetCoordim(), Nektar::StdRegions::StdExpansion::GetEdgeNormal(), Nektar::StdRegions::StdExpansion::GetEdgeToElementMap(), Nektar::StdRegions::StdExpansion::GetEorient(), Nektar::LocalRegions::Expansion::GetLocMatrix(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), GetPhysEdgeVarCoeffsFromElement(), m_negatedNormals, Vmath::Neg(), Nektar::StdRegions::NullConstFactorMap, sign, and Vmath::Vmul().

Referenced by AddHDGHelmholtzTraceTerms().

{
int i, j, n;
int nquad_e = EdgeExp[edge]->GetNumPoints(0);
int order_e = EdgeExp[edge]->GetNcoeffs();
int coordim = GetCoordim();
int ncoeffs = GetNcoeffs();
Array<OneD, NekDouble> inval (nquad_e);
Array<OneD, NekDouble> outcoeff(order_e);
Array<OneD, NekDouble> tmpcoeff(ncoeffs);
const Array<OneD, const Array<OneD, NekDouble> > &normals
= GetEdgeNormal(edge);
Array<OneD,unsigned int> emap;
Array<OneD,int> sign;
DNekVec Coeffs (ncoeffs,outarray,eWrapper);
DNekVec Tmpcoeff(ncoeffs,tmpcoeff,eWrapper);
GetEdgeToElementMap(edge,edgedir,emap,sign);
StdRegions::VarCoeffMap::const_iterator x;
/// @TODO: What direction to use here??
if ((x = varcoeffs.find(VarCoeff[0])) != varcoeffs.end())
{
Array<OneD, NekDouble> work(nquad_e);
edge, EdgeExp[edge], x->second, work);
Vmath::Vmul(nquad_e, work, 1, edgePhys, 1, edgePhys, 1);
}
//================================================================
// Add F = \tau <phi_i,in_phys>
// Fill edge and take inner product
EdgeExp[edge]->IProductWRTBase(edgePhys, outcoeff);
// add data to out array
for(i = 0; i < order_e; ++i)
{
outarray[emap[i]] += sign[i] * tau * outcoeff[i];
}
//================================================================
//===============================================================
// Add -\sum_i D_i^T M^{-1} G_i + E_i M^{-1} G_i =
// \sum_i D_i M^{-1} G_i term
// Two independent direction
for(n = 0; n < coordim; ++n)
{
Vmath::Vmul(nquad_e, normals[n], 1, edgePhys, 1, inval, 1);
if (m_negatedNormals[edge])
{
Vmath::Neg(nquad_e, inval, 1);
}
// Multiply by variable coefficient
/// @TODO: Document this (probably not needed)
// StdRegions::VarCoeffMap::const_iterator x;
// if ((x = varcoeffs.find(VarCoeff[n])) != varcoeffs.end())
// {
// GetPhysEdgeVarCoeffsFromElement(edge,EdgeExp[edge],x->second,varcoeff_work);
// Vmath::Vmul(nquad_e,varcoeff_work,1,EdgeExp[edge]->GetPhys(),1,EdgeExp[edge]->UpdatePhys(),1);
// }
EdgeExp[edge]->IProductWRTBase(inval, outcoeff);
// M^{-1} G
for(i = 0; i < ncoeffs; ++i)
{
tmpcoeff[i] = 0;
for(j = 0; j < order_e; ++j)
{
tmpcoeff[i] += invMass(i,emap[j])*sign[j]*outcoeff[j];
}
}
if(varcoeffs.find(VarCoeff[n]) != varcoeffs.end())
{
MatrixKey mkey(DerivType[n], DetShapeType(), *this, StdRegions::NullConstFactorMap, varcoeffs);
DNekScalMat &Dmat = *GetLocMatrix(mkey);
Coeffs = Coeffs + Dmat*Tmpcoeff;
}
else
{
DNekScalMat &Dmat = *GetLocMatrix(DerivType[n]);
Coeffs = Coeffs + Dmat*Tmpcoeff;
}
}
}
void Nektar::LocalRegions::Expansion2D::AddHDGHelmholtzTraceTerms ( const NekDouble  tau,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, StdRegions::StdExpansionSharedPtr > &  EdgeExp,
const StdRegions::VarCoeffMap dirForcing,
Array< OneD, NekDouble > &  outarray 
)
inline

Definition at line 459 of file Expansion2D.cpp.

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

Referenced by v_GenMatrix().

{
ASSERTL0(&inarray[0] != &outarray[0],
"Input and output arrays use the same memory");
int e, cnt, order_e, nedges = GetNedges();
Array<OneD, const NekDouble> tmp;
cnt = 0;
for(e = 0; e < nedges; ++e)
{
order_e = EdgeExp[e]->GetNcoeffs();
Array<OneD, NekDouble> edgeCoeffs(order_e);
Array<OneD, NekDouble> edgePhys (EdgeExp[e]->GetTotPoints());
Vmath::Vcopy(order_e, tmp = inarray + cnt, 1, edgeCoeffs, 1);
EdgeExp[e]->BwdTrans(edgeCoeffs, edgePhys);
tau, e, EdgeExp, edgePhys, dirForcing, outarray);
cnt += order_e;
}
}
void Nektar::LocalRegions::Expansion2D::AddNormTraceInt ( const int  dir,
Array< OneD, StdRegions::StdExpansionSharedPtr > &  EdgeExp,
Array< OneD, Array< OneD, NekDouble > > &  edgeCoeffs,
Array< OneD, NekDouble > &  outarray 
)
inline

Definition at line 382 of file Expansion2D.cpp.

References AddEdgeBoundaryInt(), Nektar::StdRegions::StdExpansion::GetEdgeNormal(), Nektar::StdRegions::StdExpansion::GetNedges(), m_negatedNormals, Vmath::Neg(), and Vmath::Vmul().

Referenced by v_DGDeriv(), and v_GenMatrix().

{
int e;
int nquad_e;
int nedges = GetNedges();
for(e = 0; e < nedges; ++e)
{
nquad_e = EdgeExp[e]->GetNumPoints(0);
Array<OneD, NekDouble> edgePhys(nquad_e);
const Array<OneD, const Array<OneD, NekDouble> > &normals
EdgeExp[e]->BwdTrans(edgeCoeffs[e], edgePhys);
Vmath::Vmul(nquad_e, normals[dir], 1, edgePhys, 1, edgePhys, 1);
{
Vmath::Neg(nquad_e, edgePhys, 1);
}
AddEdgeBoundaryInt(e, EdgeExp[e], edgePhys, outarray);
}
}
void Nektar::LocalRegions::Expansion2D::AddNormTraceInt ( const int  dir,
Array< OneD, const NekDouble > &  inarray,
Array< OneD, StdRegions::StdExpansionSharedPtr > &  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.

: Document this

Definition at line 327 of file Expansion2D.cpp.

References AddEdgeBoundaryInt(), Nektar::StdRegions::StdExpansion::GetEdgeNormal(), Nektar::StdRegions::StdExpansion::GetNedges(), m_negatedNormals, Vmath::Neg(), and Vmath::Vmul().

{
int i,e,cnt;
int order_e,nquad_e;
int nedges = GetNedges();
cnt = 0;
for(e = 0; e < nedges; ++e)
{
order_e = EdgeExp[e]->GetNcoeffs();
nquad_e = EdgeExp[e]->GetNumPoints(0);
const Array<OneD, const Array<OneD, NekDouble> > &normals
Array<OneD, NekDouble> edgeCoeffs(order_e);
Array<OneD, NekDouble> edgePhys (nquad_e);
for(i = 0; i < order_e; ++i)
{
edgeCoeffs[i] = inarray[i+cnt];
}
cnt += order_e;
EdgeExp[e]->BwdTrans(edgeCoeffs, edgePhys);
// Multiply by variable coefficient
/// @TODO: Document this
// StdRegions::VarCoeffType VarCoeff[3] = {StdRegions::eVarCoeffD00,
// StdRegions::eVarCoeffD11,
// StdRegions::eVarCoeffD22};
// StdRegions::VarCoeffMap::const_iterator x;
// Array<OneD, NekDouble> varcoeff_work(nquad_e);
// if ((x = varcoeffs.find(VarCoeff[dir])) != varcoeffs.end())
// {
// GetPhysEdgeVarCoeffsFromElement(e,EdgeExp[e],x->second,varcoeff_work);
// Vmath::Vmul(nquad_e,varcoeff_work,1,EdgeExp[e]->GetPhys(),1,EdgeExp[e]->UpdatePhys(),1);
// }
Vmath::Vmul(nquad_e, normals[dir], 1, edgePhys, 1, edgePhys, 1);
{
Vmath::Neg(nquad_e, edgePhys, 1);
}
AddEdgeBoundaryInt(e, EdgeExp[e], edgePhys, outarray, varcoeffs);
}
}
ExpansionSharedPtr Nektar::LocalRegions::Expansion2D::GetEdgeExp ( int  edge,
bool  SetUpNormal = true 
)
inline

Definition at line 175 of file Expansion2D.h.

References ASSERTL1, Nektar::StdRegions::StdExpansion::GetNedges(), and m_edgeExp.

Referenced by v_GenMatrix().

{
ASSERTL1(edge < GetNedges(), "Edge out of range.");
return m_edgeExp[edge].lock();
}
SpatialDomains::Geometry2DSharedPtr Nektar::LocalRegions::Expansion2D::GetGeom2D ( ) const
inline

Definition at line 230 of file Expansion2D.h.

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

Referenced by Nektar::MultiRegions::ExpList1D::ExpList1D(), Nektar::MultiRegions::ExpList2D::ExpList2D(), Nektar::MultiRegions::AssemblyMapCG3D::SetUp3DExpansionC0ContMap(), Nektar::MultiRegions::AssemblyMapDG::SetUpUniversalDGMap(), v_AddEdgeNormBoundaryInt(), Nektar::LocalRegions::TriExp::v_FwdTrans_BndConstrained(), Nektar::LocalRegions::QuadExp::v_FwdTrans_BndConstrained(), Nektar::LocalRegions::TriExp::v_GetCartesianEorient(), Nektar::LocalRegions::QuadExp::v_GetCartesianEorient(), v_GetEdgeInverseBoundaryMap(), and Nektar::LocalRegions::TriExp::v_GetEorient().

{
return boost::dynamic_pointer_cast<SpatialDomains::Geometry2D>(m_geom);
}
Expansion3DSharedPtr Nektar::LocalRegions::Expansion2D::GetLeftAdjacentElementExp ( ) const
inline

Definition at line 192 of file Expansion2D.h.

References ASSERTL1, and m_elementLeft.

Referenced by Nektar::MultiRegions::ExpList2D::v_GetNormals(), Nektar::LocalRegions::TriExp::v_NormVectorIProductWRTBase(), and Nektar::LocalRegions::QuadExp::v_NormVectorIProductWRTBase().

{
ASSERTL1(m_elementLeft.lock().get(), "Left adjacent element not set.");
return m_elementLeft.lock();
}
int Nektar::LocalRegions::Expansion2D::GetLeftAdjacentElementFace ( ) const
inline

Definition at line 204 of file Expansion2D.h.

References m_elementFaceLeft.

Referenced by Nektar::LocalRegions::TriExp::v_NormVectorIProductWRTBase(), and Nektar::LocalRegions::QuadExp::v_NormVectorIProductWRTBase().

{
}
void Nektar::LocalRegions::Expansion2D::GetPhysEdgeVarCoeffsFromElement ( const int  edge,
StdRegions::StdExpansionSharedPtr EdgeExp,
const Array< OneD, const NekDouble > &  varcoeff,
Array< OneD, NekDouble > &  outarray 
)
protected

Extracts the variable coefficients along an edge

Definition at line 606 of file Expansion2D.cpp.

References Nektar::StdRegions::StdExpansion::FwdTrans(), Nektar::StdRegions::StdExpansion::GetEdgeToElementMap(), Nektar::StdRegions::StdExpansion::GetEorient(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), and sign.

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

{
Array<OneD, NekDouble> tmp(GetNcoeffs());
Array<OneD, NekDouble> edgetmp(EdgeExp->GetNcoeffs());
// FwdTrans varcoeffs
FwdTrans(varcoeff, tmp);
// Map to edge
Array<OneD,unsigned int> emap;
Array<OneD, int> sign;
GetEdgeToElementMap(edge,edgedir,emap,sign);
for (unsigned int i = 0; i < EdgeExp->GetNcoeffs(); ++i)
{
edgetmp[i] = tmp[emap[i]];
}
// BwdTrans
EdgeExp->BwdTrans(edgetmp, outarray);
}
Expansion3DSharedPtr Nektar::LocalRegions::Expansion2D::GetRightAdjacentElementExp ( ) const
inline

Definition at line 198 of file Expansion2D.h.

References ASSERTL1, m_elementLeft, and m_elementRight.

{
ASSERTL1(m_elementLeft.lock().get(), "Right adjacent element not set.");
return m_elementRight.lock();
}
int Nektar::LocalRegions::Expansion2D::GetRightAdjacentElementFace ( ) const
inline

Definition at line 209 of file Expansion2D.h.

References m_elementFaceRight.

{
}
void Nektar::LocalRegions::Expansion2D::SetAdjacentElementExp ( int  face,
Expansion3DSharedPtr f 
)
inline

Definition at line 214 of file Expansion2D.h.

References ASSERTL1, m_elementFaceLeft, m_elementFaceRight, m_elementLeft, and m_elementRight.

Referenced by Nektar::MultiRegions::DisContField3D::SetUpDG().

{
if (m_elementLeft.lock().get())
{
ASSERTL1(!m_elementRight.lock().get(),
"Both adjacent elements already set.");
}
else
{
}
}
void Nektar::LocalRegions::Expansion2D::SetEdgeExp ( const int  edge,
ExpansionSharedPtr e 
)
inline

Definition at line 181 of file Expansion2D.h.

References ASSERTL1, Nektar::StdRegions::StdExpansion::GetNedges(), and m_edgeExp.

{
unsigned int nEdges = GetNedges();
ASSERTL1(edge < nEdges, "Edge out of range.");
if (m_edgeExp.size() < nEdges)
{
m_edgeExp.resize(nEdges);
}
m_edgeExp[edge] = e;
}
void Nektar::LocalRegions::Expansion2D::SetTraceToGeomOrientation ( Array< OneD, StdRegions::StdExpansionSharedPtr > &  EdgeExp,
Array< OneD, NekDouble > &  inout 
)

Definition at line 306 of file Expansion2D.cpp.

References Nektar::StdRegions::StdExpansion::GetEdgeNcoeffs(), Nektar::StdRegions::StdExpansion::GetEorient(), and Nektar::StdRegions::StdExpansion::GetNedges().

Referenced by v_GenMatrix().

{
int i, cnt = 0;
int nedges = GetNedges();
Array<OneD, NekDouble> e_tmp;
for(i = 0; i < nedges; ++i)
{
EdgeExp[i]->SetCoeffsToOrientation(GetEorient(i),
e_tmp = inout + cnt,
e_tmp = inout + cnt);
cnt += GetEdgeNcoeffs(i);
}
}
void Nektar::LocalRegions::Expansion2D::v_AddEdgeNormBoundaryInt ( const int  edge,
StdRegions::StdExpansionSharedPtr EdgeExp,
const Array< OneD, const NekDouble > &  Fx,
const Array< OneD, const NekDouble > &  Fy,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Definition at line 54 of file Expansion2D.cpp.

References Nektar::StdRegions::StdExpansion::AddEdgeNormBoundaryInt(), ASSERTL1, Nektar::StdRegions::StdExpansion::GetCoordim(), Nektar::StdRegions::StdExpansion::GetEdgeNormal(), GetGeom2D(), Nektar::StdRegions::StdExpansion::GetNedges(), m_edgeExp, m_negatedNormals, m_requireNeg, Vmath::Neg(), Vmath::Vmul(), and Vmath::Vvtvp().

{
"Routine only set up for two-dimensions");
const Array<OneD, const Array<OneD, NekDouble> > normals
= GetEdgeNormal(edge);
if (m_requireNeg.size() == 0)
{
for (int i = 0; i < GetNedges(); ++i)
{
m_requireNeg[i] = false;
{
m_requireNeg[i] = true;
continue;
}
m_edgeExp[i].lock()->as<Expansion1D>();
if (edgeExp->GetRightAdjacentElementExp())
{
if (edgeExp->GetRightAdjacentElementExp()->GetGeom2D()
->GetGlobalID() == GetGeom2D()->GetGlobalID())
{
m_requireNeg[i] = true;
}
}
}
}
// We allow the case of mixed polynomial order by supporting only
// those modes on the edge common to both adjoining elements. This
// is enforced here by taking the minimum size and padding with
// zeros.
int nquad_e = min(EdgeExp->GetNumPoints(0),
int(normals[0].num_elements()));
int nEdgePts = EdgeExp->GetTotPoints();
Array<OneD, NekDouble> edgePhys(nEdgePts);
Vmath::Vmul (nquad_e, normals[0], 1, Fx, 1, edgePhys, 1);
Vmath::Vvtvp(nquad_e, normals[1], 1, Fy, 1, edgePhys, 1,
edgePhys, 1);
Expansion1DSharedPtr locExp = EdgeExp->as<Expansion1D>();
if (m_negatedNormals[edge])
{
Vmath::Neg(nquad_e, edgePhys, 1);
}
else if (locExp->GetRightAdjacentElementEdge() != -1)
{
if (locExp->GetRightAdjacentElementExp()->GetGeom2D()->GetGlobalID()
== GetGeom2D()->GetGlobalID())
{
Vmath::Neg(nquad_e, edgePhys, 1);
}
}
AddEdgeNormBoundaryInt(edge, EdgeExp, edgePhys, outarray);
}
void Nektar::LocalRegions::Expansion2D::v_AddEdgeNormBoundaryInt ( const int  edge,
StdRegions::StdExpansionSharedPtr EdgeExp,
const Array< OneD, const NekDouble > &  Fn,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Definition at line 127 of file Expansion2D.cpp.

References ASSERTL0, Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::StdRegions::eEdgeToElement, Nektar::StdRegions::eFactorGaussEdge, Nektar::LibUtilities::eGauss_Lagrange, Nektar::StdRegions::eGaussDG, Nektar::StdRegions::eMass, Nektar::LibUtilities::eSegment, Nektar::StdRegions::StdExpansion::GetBasisNumModes(), Nektar::StdRegions::StdExpansion::GetEorient(), GetGeom2D(), Nektar::StdRegions::StdExpansion::GetIndexMap(), Nektar::StdRegions::StdExpansion::GetNedges(), Nektar::LibUtilities::InterpCoeff1D(), Nektar::StdRegions::StdExpansion::m_base, m_edgeExp, m_negatedNormals, m_requireNeg, Nektar::StdRegions::StdExpansion::m_stdMatrixManager, and Vmath::Neg().

{
int i;
if (m_requireNeg.size() == 0)
{
for (i = 0; i < GetNedges(); ++i)
{
m_requireNeg[i] = false;
{
m_requireNeg[i] = true;
continue;
}
m_edgeExp[i].lock()->as<Expansion1D>();
if (edgeExp->GetRightAdjacentElementExp())
{
if (edgeExp->GetRightAdjacentElementExp()->GetGeom2D()
->GetGlobalID() == GetGeom2D()->GetGlobalID())
{
m_requireNeg[i] = true;
}
}
}
}
StdRegions::IndexMapKey ikey(
edge, GetEorient(edge));
// Order of the element
int order_e = map->num_elements();
// Order of the trace
int n_coeffs = EdgeExp->GetNcoeffs();
Array<OneD, NekDouble> edgeCoeffs(n_coeffs);
if(n_coeffs!=order_e) // Going to orthogonal space
{
EdgeExp->FwdTrans(Fn, edgeCoeffs);
Expansion1DSharedPtr locExp = EdgeExp->as<Expansion1D>();
if (m_requireNeg[edge])
{
Vmath::Neg(n_coeffs, edgeCoeffs, 1);
}
Array<OneD, NekDouble> coeff(n_coeffs,0.0);
LibUtilities::BasisType btype = ((LibUtilities::BasisType) 1); //1-->Ortho_A
LibUtilities::BasisKey bkey_ortho(btype,EdgeExp->GetBasis(0)->GetNumModes(),EdgeExp->GetBasis(0)->GetPointsKey());
LibUtilities::BasisKey bkey(EdgeExp->GetBasis(0)->GetBasisType(),EdgeExp->GetBasis(0)->GetNumModes(),EdgeExp->GetBasis(0)->GetPointsKey());
LibUtilities::InterpCoeff1D(bkey,edgeCoeffs,bkey_ortho,coeff);
// Cutting high frequencies
for(i = order_e; i < n_coeffs; i++)
{
coeff[i] = 0.0;
}
LibUtilities::InterpCoeff1D(bkey_ortho,coeff,bkey,edgeCoeffs);
StdRegions::StdMatrixKey masskey(StdRegions::eMass,LibUtilities::eSegment,*EdgeExp);
EdgeExp->MassMatrixOp(edgeCoeffs, edgeCoeffs, masskey);
}
else
{
EdgeExp->IProductWRTBase(Fn, edgeCoeffs);
Expansion1DSharedPtr locExp = EdgeExp->as<Expansion1D>();
if (m_requireNeg[edge])
{
Vmath::Neg(n_coeffs, edgeCoeffs, 1);
}
}
// Implementation for all the basis except Gauss points
if(EdgeExp->GetBasis(0)->GetBasisType() !=
{
// add data to outarray if forward edge normal is outwards
for(i = 0; i < order_e; ++i)
{
outarray[(*map)[i].index] +=
(*map)[i].sign * edgeCoeffs[i];
}
}
else
{
int nCoeffs0, nCoeffs1;
int j;
StdRegions::StdMatrixKey key(StdRegions::eGaussDG,
DetShapeType(),*this,factors);
switch(edge)
{
case 0:
{
nCoeffs1 = m_base[1]->GetNumModes();
for(i = 0; i < order_e; ++i)
{
for(j = 0; j < nCoeffs1; j++)
{
outarray[(*map)[i].index + j*order_e] +=
mat_gauss->GetPtr()[j]*
(*map)[i].sign*edgeCoeffs[i];
}
}
break;
}
case 1:
{
nCoeffs0 = m_base[0]->GetNumModes();
for(i = 0; i < order_e; ++i)
{
for(j = 0; j < nCoeffs0; j++)
{
outarray[(*map)[i].index - j] +=
mat_gauss->GetPtr()[order_e - 1 -j]*
(*map)[i].sign*edgeCoeffs[i];
}
}
break;
}
case 2:
{
nCoeffs1 = m_base[1]->GetNumModes();
for(i = 0; i < order_e; ++i)
{
for(j = 0; j < nCoeffs1; j++)
{
outarray[(*map)[i].index - j*order_e] +=
mat_gauss->GetPtr()[order_e - 1 - j]*
(*map)[i].sign*edgeCoeffs[i];
}
}
break;
}
case 3:
{
nCoeffs0 = m_base[0]->GetNumModes();
for(i = 0; i < order_e; ++i)
{
for(j = 0; j < nCoeffs0; j++)
{
outarray[(*map)[i].index + j] +=
mat_gauss->GetPtr()[j]*
(*map)[i].sign*edgeCoeffs[i];
}
}
break;
}
default:
ASSERTL0(false,"edge value (< 3) is out of range");
break;
}
}
}
void Nektar::LocalRegions::Expansion2D::v_AddRobinEdgeContribution ( const int  edgeid,
const Array< OneD, const NekDouble > &  primCoeffs,
Array< OneD, NekDouble > &  coeffs 
)
protectedvirtual

Given an edge and vector of element coefficients:

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1260 of file Expansion2D.cpp.

References ASSERTL1, Nektar::StdRegions::eMass, Nektar::LibUtilities::eSegment, Nektar::StdRegions::eVarCoeffMass, Nektar::StdRegions::StdExpansion::GetEdgeToElementMap(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), Nektar::StdRegions::StdExpansion::IsBoundaryInteriorExpansion(), m_edgeExp, Nektar::StdRegions::NullConstFactorMap, sign, Nektar::StdRegions::StdExpansion::v_GetEorient(), and Vmath::Zero().

{
"Not set up for non boundary-interior expansions");
int i;
ExpansionSharedPtr edgeExp = m_edgeExp[edgeid].lock();
int order_e = edgeExp->GetNcoeffs();
Array<OneD,unsigned int> map;
Array<OneD,int> sign;
varcoeffs[StdRegions::eVarCoeffMass] = primCoeffs;
LocalRegions::MatrixKey mkey(StdRegions::eMass,LibUtilities::eSegment, *edgeExp, StdRegions::NullConstFactorMap, varcoeffs);
DNekScalMat &edgemat = *edgeExp->GetLocMatrix(mkey);
NekVector<NekDouble> vEdgeCoeffs (order_e);
GetEdgeToElementMap(edgeid,v_GetEorient(edgeid),map,sign);
for (i = 0; i < order_e; ++i)
{
vEdgeCoeffs[i] = coeffs[map[i]]*sign[i];
}
Vmath::Zero(GetNcoeffs(), coeffs, 1);
vEdgeCoeffs = edgemat * vEdgeCoeffs;
for (i = 0; i < order_e; ++i)
{
coeffs[map[i]] = vEdgeCoeffs[i]*sign[i];
}
}
void Nektar::LocalRegions::Expansion2D::v_AddRobinMassMatrix ( const int  edgeid,
const Array< OneD, const NekDouble > &  primCoeffs,
DNekMatSharedPtr inoutmat 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1135 of file Expansion2D.cpp.

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::GetEdgeNcoeffs(), Nektar::StdRegions::StdExpansion::GetEdgeToElementMap(), Nektar::StdRegions::StdExpansion::GetEorient(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), Nektar::StdRegions::StdExpansion::IsBoundaryInteriorExpansion(), m_edgeExp, Nektar::StdRegions::NullConstFactorMap, Nektar::StdRegions::StdExpansion::NumBndryCoeffs(), Nektar::StdRegions::StdExpansion::NumDGBndryCoeffs(), sign, and Nektar::StdRegions::StdExpansion::v_GetEorient().

{
"Not set up for non boundary-interior expansions");
ASSERTL1(inoutmat->GetRows() == inoutmat->GetColumns(),
"Assuming that input matrix was square");
int i,j;
int id1,id2;
ExpansionSharedPtr edgeExp = m_edgeExp[edge].lock();
int order_e = edgeExp->GetNcoeffs();
Array<OneD,unsigned int> map;
Array<OneD,int> sign;
varcoeffs[StdRegions::eVarCoeffMass] = primCoeffs;
LocalRegions::MatrixKey mkey(StdRegions::eMass,LibUtilities::eSegment, *edgeExp, StdRegions::NullConstFactorMap, varcoeffs);
DNekScalMat &edgemat = *edgeExp->GetLocMatrix(mkey);
// Now need to identify a map which takes the local edge
// mass matrix to the matrix stored in inoutmat;
// This can currently be deduced from the size of the matrix
// - if inoutmat.m_rows() == v_NCoeffs() it is a full
// matrix system
// - if inoutmat.m_rows() == v_NumBndCoeffs() it is a
// boundary CG system
// - if inoutmat.m_rows() == v_NumDGBndCoeffs() it is a
// trace DG system
int rows = inoutmat->GetRows();
if (rows == GetNcoeffs())
{
GetEdgeToElementMap(edge,v_GetEorient(edge),map,sign);
}
else if(rows == NumBndryCoeffs())
{
int nbndry = NumBndryCoeffs();
Array<OneD,unsigned int> bmap(nbndry);
GetEdgeToElementMap(edge,v_GetEorient(edge),map,sign);
for(i = 0; i < order_e; ++i)
{
for(j = 0; j < nbndry; ++j)
{
if(map[i] == bmap[j])
{
map[i] = j;
break;
}
}
ASSERTL1(j != nbndry,"Did not find number in map");
}
}
else if (rows == NumDGBndryCoeffs())
{
// possibly this should be a separate method
int cnt = 0;
map = Array<OneD, unsigned int> (order_e);
sign = Array<OneD, int> (order_e,1);
for(i = 0; i < edge; ++i)
{
cnt += GetEdgeNcoeffs(i);
}
for(i = 0; i < order_e; ++i)
{
map[i] = cnt++;
}
// check for mapping reversal
{
switch(edgeExp->GetBasis(0)->GetBasisType())
{
reverse( map.get() , map.get()+order_e);
break;
reverse( map.get() , map.get()+order_e);
break;
{
swap(map[0],map[1]);
for(i = 3; i < order_e; i+=2)
{
sign[i] = -1;
}
}
break;
default:
ASSERTL0(false,"Edge boundary type not valid for this method");
}
}
}
else
{
ASSERTL0(false,"Could not identify matrix type from dimension");
}
for(i = 0; i < order_e; ++i)
{
id1 = map[i];
for(j = 0; j < order_e; ++j)
{
id2 = map[j];
(*inoutmat)(id1,id2) += edgemat(i,j)*sign[i]*sign[j];
}
}
}
DNekMatSharedPtr Nektar::LocalRegions::Expansion2D::v_BuildVertexMatrix ( const DNekScalMatSharedPtr r_bnd)
protectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1295 of file Expansion2D.cpp.

References Nektar::eFULL, Nektar::StdRegions::StdExpansion::GetNverts(), and Nektar::StdRegions::StdExpansion::GetVertexMap().

{
MatrixStorage storage = eFULL;
DNekMatSharedPtr m_vertexmatrix;
int nVerts, vid1, vid2, vMap1, vMap2;
NekDouble VertexValue;
nVerts = GetNverts();
m_vertexmatrix =
MemoryManager<DNekMat>::AllocateSharedPtr(
nVerts, nVerts, 0.0, storage);
DNekMat &VertexMat = (*m_vertexmatrix);
for (vid1 = 0; vid1 < nVerts; ++vid1)
{
vMap1 = GetVertexMap(vid1);
for (vid2 = 0; vid2 < nVerts; ++vid2)
{
vMap2 = GetVertexMap(vid2);
VertexValue = (*r_bnd)(vMap1, vMap2);
VertexMat.SetValue(vid1, vid2, VertexValue);
}
}
return m_vertexmatrix;
}
void Nektar::LocalRegions::Expansion2D::v_DGDeriv ( const int  dir,
const Array< OneD, const NekDouble > &  incoeffs,
Array< OneD, StdRegions::StdExpansionSharedPtr > &  EdgeExp,
Array< OneD, Array< OneD, NekDouble > > &  edgeCoeffs,
Array< OneD, NekDouble > &  out_d 
)
protectedvirtual

Definition at line 1097 of file Expansion2D.cpp.

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

{
int ncoeffs = GetNcoeffs();
DNekScalMat &Dmat = *GetLocMatrix(DerivType[dir]);
Array<OneD, NekDouble> coeffs = incoeffs;
DNekVec Coeffs (ncoeffs,coeffs, eWrapper);
Coeffs = Transpose(Dmat)*Coeffs;
Vmath::Neg(ncoeffs, coeffs,1);
// Add the boundary integral including the relevant part of
// the normal
AddNormTraceInt(dir, EdgeExp, edgeCoeffs, coeffs);
DNekVec Out_d (ncoeffs,out_d,eWrapper);
Out_d = InvMass*Coeffs;
}
bool Nektar::LocalRegions::Expansion2D::v_EdgeNormalNegated ( const int  edge)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1395 of file Expansion2D.cpp.

References m_negatedNormals.

{
return m_negatedNormals[edge];
}
DNekMatSharedPtr Nektar::LocalRegions::Expansion2D::v_GenMatrix ( const StdRegions::StdMatrixKey mkey)
protectedvirtual

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::QuadExp, Nektar::LocalRegions::TriExp, and Nektar::LocalRegions::NodalTriExp.

Definition at line 641 of file Expansion2D.cpp.

References AddHDGHelmholtzTraceTerms(), AddNormTraceInt(), 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::eVarCoeffD00, Nektar::StdRegions::eVarCoeffD11, Nektar::StdRegions::eVarCoeffD22, Nektar::StdRegions::eVarCoeffMass, Nektar::StdRegions::eWeakDeriv0, Nektar::StdRegions::eWeakDeriv1, Nektar::StdRegions::eWeakDeriv2, Nektar::eWrapper, Vmath::Fill(), Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdMatrixKey::GetConstFactors(), Nektar::StdRegions::StdExpansion::GetCoordim(), GetEdgeExp(), Nektar::StdRegions::StdExpansion::GetEdgeNormal(), Nektar::StdRegions::StdExpansion::GetEdgeToElementMap(), Nektar::StdRegions::StdExpansion::GetEorient(), Nektar::LocalRegions::Expansion::GetLocMatrix(), Nektar::StdRegions::StdMatrixKey::GetMatrixType(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), Nektar::StdRegions::StdExpansion::GetNedges(), GetPhysEdgeVarCoeffsFromElement(), Nektar::StdRegions::StdExpansion::GetTotPoints(), Nektar::StdRegions::StdMatrixKey::GetVarCoeff(), Nektar::StdRegions::StdMatrixKey::GetVarCoeffAsMap(), Nektar::StdRegions::StdMatrixKey::GetVarCoeffs(), Nektar::StdRegions::StdMatrixKey::HasVarCoeff(), Nektar::StdRegions::StdExpansion::IProductWRTBase(), Nektar::StdRegions::StdExpansion::IsBoundaryInteriorExpansion(), Nektar::StdRegions::StdExpansion::m_ncoeffs, m_negatedNormals, Vmath::Neg(), Nektar::StdRegions::NullConstFactorMap, Nektar::StdRegions::StdExpansion::NumDGBndryCoeffs(), SetTraceToGeomOrientation(), sign, Vmath::Svtvp(), Nektar::Transpose(), Vmath::Vcopy(), Vmath::Vmul(), Vmath::Vvtvp(), and Vmath::Zero().

{
DNekMatSharedPtr returnval;
switch(mkey.GetMatrixType())
{
// (Z^e)^{-1} (Eqn. 33, P22)
{
"HybridDGHelmholtz matrix not set up "
"for non boundary-interior expansions");
int i,j,k;
NekDouble lambdaval = mkey.GetConstFactor(StdRegions::eFactorLambda);
NekDouble tau = mkey.GetConstFactor(StdRegions::eFactorTau);
int ncoeffs = GetNcoeffs();
int nedges = GetNedges();
Array<OneD,unsigned int> emap;
Array<OneD,int> sign;
int order_e, coordim = GetCoordim();
DNekMat LocMat(ncoeffs,ncoeffs);
returnval = MemoryManager<DNekMat>::AllocateSharedPtr(ncoeffs,ncoeffs);
DNekMat &Mat = *returnval;
Vmath::Zero(ncoeffs*ncoeffs,Mat.GetPtr(),1);
for(i=0; i < coordim; ++i)
{
if(mkey.HasVarCoeff(Coeffs[i]))
{
MatrixKey DmatkeyL(DerivType[i], DetShapeType(), *this, StdRegions::NullConstFactorMap, mkey.GetVarCoeffAsMap(Coeffs[i]));
MatrixKey DmatkeyR(DerivType[i], DetShapeType(), *this);
DNekScalMat &DmatL = *GetLocMatrix(DmatkeyL);
DNekScalMat &DmatR = *GetLocMatrix(DmatkeyR);
Mat = Mat + DmatL*invMass*Transpose(DmatR);
}
else
{
DNekScalMat &Dmat = *GetLocMatrix(DerivType[i]);
Mat = Mat + Dmat*invMass*Transpose(Dmat);
}
}
// Add Mass Matrix Contribution for Helmholtz problem
Mat = Mat + lambdaval*Mass;
// Add tau*E_l using elemental mass matrices on each edge
for(i = 0; i < nedges; ++i)
{
EdgeExp = GetEdgeExp(i);
EdgeExp2 = GetEdgeExp(i);
order_e = EdgeExp->GetNcoeffs();
int nq = EdgeExp->GetNumPoints(0);
GetEdgeToElementMap(i,edgedir,emap,sign);
// @TODO: Document
StdRegions::VarCoeffMap edgeVarCoeffs;
if (mkey.HasVarCoeff(StdRegions::eVarCoeffD00))
{
Array<OneD, NekDouble> mu(nq);
GetPhysEdgeVarCoeffsFromElement(i, EdgeExp2, mkey.GetVarCoeff(StdRegions::eVarCoeffD00), mu);
edgeVarCoeffs[StdRegions::eVarCoeffMass] = mu;
}
DNekScalMat &eMass = *EdgeExp->GetLocMatrix(StdRegions::eMass, StdRegions::NullConstFactorMap, edgeVarCoeffs);
//DNekScalMat &eMass = *EdgeExp->GetLocMatrix(StdRegions::eMass);
for(j = 0; j < order_e; ++j)
{
for(k = 0; k < order_e; ++k)
{
Mat(emap[j],emap[k]) = Mat(emap[j],emap[k]) + tau*sign[j]*sign[k]*eMass(j,k);
}
}
}
}
break;
// U^e (P22)
{
int i,j,k;
int nbndry = NumDGBndryCoeffs();
int ncoeffs = GetNcoeffs();
int nedges = GetNedges();
NekDouble tau = mkey.GetConstFactor(StdRegions::eFactorTau);
Array<OneD,NekDouble> lambda(nbndry);
DNekVec Lambda(nbndry,lambda,eWrapper);
Array<OneD,NekDouble> ulam(ncoeffs);
DNekVec Ulam(ncoeffs,ulam,eWrapper);
Array<OneD,NekDouble> f(ncoeffs);
DNekVec F(ncoeffs,f,eWrapper);
Array<OneD,StdRegions::StdExpansionSharedPtr> EdgeExp(nedges);
// declare matrix space
returnval = MemoryManager<DNekMat>::AllocateSharedPtr(ncoeffs,nbndry);
DNekMat &Umat = *returnval;
// Z^e matrix
MatrixKey newkey(StdRegions::eInvHybridDGHelmholtz, DetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
DNekScalMat &invHmat = *GetLocMatrix(newkey);
Array<OneD,unsigned int> emap;
Array<OneD,int> sign;
for(i = 0; i < nedges; ++i)
{
EdgeExp[i] = GetEdgeExp(i);
}
// for each degree of freedom of the lambda space
// calculate Umat entry
// Generate Lambda to U_lambda matrix
for(j = 0; j < nbndry; ++j)
{
// standard basis vectors e_j
Vmath::Zero(nbndry,&lambda[0],1);
Vmath::Zero(ncoeffs,&f[0],1);
lambda[j] = 1.0;
SetTraceToGeomOrientation(EdgeExp,lambda);
// Compute F = [I D_1 M^{-1} D_2 M^{-1}] C e_j
tau, lambda, EdgeExp, mkey.GetVarCoeffs(), f);
// Compute U^e_j
Ulam = invHmat*F; // generate Ulam from lambda
// fill column of matrix
for(k = 0; k < ncoeffs; ++k)
{
Umat(k,j) = Ulam[k];
}
}
}
break;
// Q_0, Q_1, Q_2 matrices (P23)
// Each are a product of a row of Eqn 32 with the C matrix.
// Rather than explicitly computing all of Eqn 32, we note each
// row is almost a multiple of U^e, so use that as our starting
// point.
{
int i,j,k,dir;
int nbndry = NumDGBndryCoeffs();
int ncoeffs = GetNcoeffs();
int nedges = GetNedges();
Array<OneD,NekDouble> lambda(nbndry);
DNekVec Lambda(nbndry,lambda,eWrapper);
Array<OneD,StdRegions::StdExpansionSharedPtr> EdgeExp(nedges);
Array<OneD,NekDouble> ulam(ncoeffs);
DNekVec Ulam(ncoeffs,ulam,eWrapper);
Array<OneD,NekDouble> f(ncoeffs);
DNekVec F(ncoeffs,f,eWrapper);
// declare matrix space
returnval = MemoryManager<DNekMat>::AllocateSharedPtr(ncoeffs,nbndry);
DNekMat &Qmat = *returnval;
// Lambda to U matrix
MatrixKey lamToUkey(StdRegions::eHybridDGLamToU, DetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
DNekScalMat &lamToU = *GetLocMatrix(lamToUkey);
// Inverse mass matrix
for(i = 0; i < nedges; ++i)
{
EdgeExp[i] = GetEdgeExp(i);
}
//Weak Derivative matrix
switch(mkey.GetMatrixType())
{
dir = 0;
break;
dir = 1;
break;
dir = 2;
break;
default:
ASSERTL0(false,"Direction not known");
break;
}
// for each degree of freedom of the lambda space
// calculate Qmat entry
// Generate Lambda to Q_lambda matrix
for(j = 0; j < nbndry; ++j)
{
Vmath::Zero(nbndry,&lambda[0],1);
lambda[j] = 1.0;
// for lambda[j] = 1 this is the solution to ulam
for(k = 0; k < ncoeffs; ++k)
{
Ulam[k] = lamToU(k,j);
}
// -D^T ulam
Vmath::Neg(ncoeffs,&ulam[0],1);
F = Transpose(*Dmat)*Ulam;
SetTraceToGeomOrientation(EdgeExp,lambda);
// Add the C terms resulting from the I's on the
// diagonals of Eqn 32
AddNormTraceInt(dir,lambda,EdgeExp,f,mkey.GetVarCoeffs());
// finally multiply by inverse mass matrix
Ulam = invMass*F;
// fill column of matrix (Qmat is in column major format)
Vmath::Vcopy(ncoeffs,&ulam[0],1,&(Qmat.GetPtr())[0]+j*ncoeffs,1);
}
}
break;
// Matrix K (P23)
{
int i,j,e,cnt;
int order_e, nquad_e;
int nbndry = NumDGBndryCoeffs();
int coordim = GetCoordim();
int nedges = GetNedges();
NekDouble tau = mkey.GetConstFactor(StdRegions::eFactorTau);
Array<OneD,NekDouble> work, varcoeff_work;
Array<OneD,const Array<OneD, NekDouble> > normals;
Array<OneD,StdRegions::StdExpansionSharedPtr> EdgeExp(nedges);
Array<OneD, NekDouble> lam(nbndry);
Array<OneD,unsigned int> emap;
Array<OneD, int> sign;
// declare matrix space
returnval = MemoryManager<DNekMat>::AllocateSharedPtr(nbndry, nbndry);
DNekMat &BndMat = *returnval;
// Matrix to map Lambda to U
MatrixKey LamToUkey(StdRegions::eHybridDGLamToU, DetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
DNekScalMat &LamToU = *GetLocMatrix(LamToUkey);
// Matrix to map Lambda to Q0
MatrixKey LamToQ0key(StdRegions::eHybridDGLamToQ0, DetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
LamToQ[0] = GetLocMatrix(LamToQ0key);
// Matrix to map Lambda to Q1
MatrixKey LamToQ1key(StdRegions::eHybridDGLamToQ1, DetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
LamToQ[1] = GetLocMatrix(LamToQ1key);
// Matrix to map Lambda to Q2 for 3D coordinates
if (coordim == 3)
{
MatrixKey LamToQ2key(StdRegions::eHybridDGLamToQ2, DetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
LamToQ[2] = GetLocMatrix(LamToQ2key);
}
// Set up edge segment expansions from local geom info
for(i = 0; i < nedges; ++i)
{
EdgeExp[i] = GetEdgeExp(i);
}
// Set up matrix derived from <mu, Q_lam.n - \tau (U_lam - Lam) >
for(i = 0; i < nbndry; ++i)
{
cnt = 0;
Vmath::Zero(nbndry,lam,1);
lam[i] = 1.0;
for(e = 0; e < nedges; ++e)
{
order_e = EdgeExp[e]->GetNcoeffs();
nquad_e = EdgeExp[e]->GetNumPoints(0);
normals = GetEdgeNormal(e);
edgedir = GetEorient(e);
work = Array<OneD,NekDouble>(nquad_e);
varcoeff_work = Array<OneD, NekDouble>(nquad_e);
GetEdgeToElementMap(e,edgedir,emap,sign);
const StdRegions::VarCoeffMap &varcoeffs = mkey.GetVarCoeffs();
StdRegions::VarCoeffMap::const_iterator x;
// Q0 * n0 (BQ_0 terms)
Array<OneD, NekDouble> edgeCoeffs(order_e);
Array<OneD, NekDouble> edgePhys (nquad_e);
for(j = 0; j < order_e; ++j)
{
edgeCoeffs[j] = sign[j]*(*LamToQ[0])(emap[j],i);
}
EdgeExp[e]->BwdTrans(edgeCoeffs, edgePhys);
// @TODO Var coeffs
// Multiply by variable coefficient
// if ((x = varcoeffs.find(VarCoeff[0])) != varcoeffs.end())
// {
// GetPhysEdgeVarCoeffsFromElement(e,EdgeExp[e],x->second,varcoeff_work);
// Vmath::Vmul(nquad_e,varcoeff_work,1,EdgeExp[e]->GetPhys(),1,EdgeExp[e]->UpdatePhys(),1);
// }
Vmath::Vmul(nquad_e, normals[0], 1, edgePhys, 1, work, 1);
// Q1 * n1 (BQ_1 terms)
for(j = 0; j < order_e; ++j)
{
edgeCoeffs[j] = sign[j]*(*LamToQ[1])(emap[j],i);
}
EdgeExp[e]->BwdTrans(edgeCoeffs, edgePhys);
// @TODO var coeffs
// Multiply by variable coefficients
// if ((x = varcoeffs.find(VarCoeff[1])) != varcoeffs.end())
// {
// GetPhysEdgeVarCoeffsFromElement(e,EdgeExp[e],x->second,varcoeff_work);
// Vmath::Vmul(nquad_e,varcoeff_work,1,EdgeExp[e]->GetPhys(),1,EdgeExp[e]->UpdatePhys(),1);
// }
Vmath::Vvtvp(nquad_e, normals[1], 1, edgePhys, 1,
work, 1, work, 1);
// Q2 * n2 (BQ_2 terms)
if (coordim == 3)
{
for(j = 0; j < order_e; ++j)
{
edgeCoeffs[j] = sign[j]*(*LamToQ[2])(emap[j],i);
}
EdgeExp[e]->BwdTrans(edgeCoeffs, edgePhys);
// @TODO var coeffs
// Multiply by variable coefficients
// if ((x = varcoeffs.find(VarCoeff[2])) != varcoeffs.end())
// {
// GetPhysEdgeVarCoeffsFromElement(e,EdgeExp[e],x->second,varcoeff_work);
// Vmath::Vmul(nquad_e,varcoeff_work,1,EdgeExp[e]->GetPhys(),1,EdgeExp[e]->UpdatePhys(),1);
// }
Vmath::Vvtvp(nquad_e, normals[2], 1, edgePhys, 1,
work, 1, work, 1);
}
{
Vmath::Neg(nquad_e, work, 1);
}
// - tau (ulam - lam)
// Corresponds to the G and BU terms.
for(j = 0; j < order_e; ++j)
{
edgeCoeffs[j] = sign[j]*LamToU(emap[j],i) - lam[cnt+j];
}
EdgeExp[e]->BwdTrans(edgeCoeffs, edgePhys);
// Multiply by variable coefficients
if ((x = varcoeffs.find(VarCoeff[0])) != varcoeffs.end())
{
GetPhysEdgeVarCoeffsFromElement(e,EdgeExp[e],x->second,varcoeff_work);
Vmath::Vmul(nquad_e,varcoeff_work,1,edgePhys,1,edgePhys,1);
}
Vmath::Svtvp(nquad_e,-tau,edgePhys,1,
work,1,work,1);
/// TODO: Add variable coeffs
EdgeExp[e]->IProductWRTBase(work, edgeCoeffs);
EdgeExp[e]->SetCoeffsToOrientation(edgeCoeffs, edgedir);
for(j = 0; j < order_e; ++j)
{
BndMat(cnt+j,i) = edgeCoeffs[j];
}
cnt += order_e;
}
}
}
break;
//HDG postprocessing
{
MatrixKey lapkey(StdRegions::eLaplacian, DetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
DNekScalMat &LapMat = *GetLocMatrix(lapkey);
returnval = MemoryManager<DNekMat>::AllocateSharedPtr(LapMat.GetRows(),LapMat.GetColumns());
DNekMatSharedPtr lmat = returnval;
(*lmat) = LapMat;
// replace first column with inner product wrt 1
int nq = GetTotPoints();
Array<OneD, NekDouble> tmp(nq);
Array<OneD, NekDouble> outarray(m_ncoeffs);
Vmath::Fill(nq,1.0,tmp,1);
IProductWRTBase(tmp, outarray);
Vmath::Vcopy(m_ncoeffs,&outarray[0],1,
&(lmat->GetPtr())[0],1);
lmat->Invert();
}
break;
default:
ASSERTL0(false,"This matrix type cannot be generated from this class");
break;
}
return returnval;
}
Array< OneD, unsigned int > Nektar::LocalRegions::Expansion2D::v_GetEdgeInverseBoundaryMap ( int  eid)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1326 of file Expansion2D.cpp.

References Nektar::StdRegions::StdExpansion::GetBoundaryMap(), Nektar::StdRegions::StdExpansion::GetEdgeInteriorMap(), Nektar::StdRegions::StdExpansion::GetEdgeNcoeffs(), GetGeom2D(), and Nektar::StdRegions::StdExpansion::NumBndryCoeffs().

{
int n, j;
int nEdgeCoeffs;
int nBndCoeffs = NumBndryCoeffs();
Array<OneD, unsigned int> bmap(nBndCoeffs);
// Map from full system to statically condensed system (i.e reverse
// GetBoundaryMap)
map<int, int> invmap;
for (j = 0; j < nBndCoeffs; ++j)
{
invmap[bmap[j]] = j;
}
// Number of interior edge coefficients
nEdgeCoeffs = GetEdgeNcoeffs(eid) - 2;
Array<OneD, unsigned int> edgemaparray(nEdgeCoeffs);
Array<OneD, unsigned int> maparray (nEdgeCoeffs);
Array<OneD, int> signarray (nEdgeCoeffs, 1);
StdRegions::Orientation eOrient = geom->GetEorient(eid);
// maparray is the location of the edge within the matrix
GetEdgeInteriorMap(eid, eOrient, maparray, signarray);
for (n = 0; n < nEdgeCoeffs; ++n)
{
edgemaparray[n] = invmap[maparray[n]];
}
return edgemaparray;
}
const StdRegions::NormalVector & Nektar::LocalRegions::Expansion2D::v_GetEdgeNormal ( const int  edge) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1370 of file Expansion2D.cpp.

References ASSERTL0, and m_edgeNormals.

Referenced by v_GetSurfaceNormal().

{
std::map<int, StdRegions::NormalVector>::const_iterator x;
x = m_edgeNormals.find(edge);
ASSERTL0 (x != m_edgeNormals.end(),
"Edge normal not computed.");
return x->second;
}
const StdRegions::NormalVector & Nektar::LocalRegions::Expansion2D::v_GetSurfaceNormal ( const int  id) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1379 of file Expansion2D.cpp.

References v_GetEdgeNormal().

{
return v_GetEdgeNormal(id);
}
void Nektar::LocalRegions::Expansion2D::v_NegateEdgeNormal ( const int  edge)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1385 of file Expansion2D.cpp.

References Nektar::StdRegions::StdExpansion::GetCoordim(), m_edgeNormals, m_negatedNormals, and Vmath::Neg().

{
m_negatedNormals[edge] = true;
for (int i = 0; i < GetCoordim(); ++i)
{
Vmath::Neg(m_edgeNormals[edge][i].num_elements(),
m_edgeNormals[edge][i], 1);
}
}
void Nektar::LocalRegions::Expansion2D::v_SetUpPhysNormals ( const int  edge)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1365 of file Expansion2D.cpp.

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

{
}

Member Data Documentation

std::vector<ExpansionWeakPtr> Nektar::LocalRegions::Expansion2D::m_edgeExp
protected

Definition at line 121 of file Expansion2D.h.

Referenced by GetEdgeExp(), SetEdgeExp(), v_AddEdgeNormBoundaryInt(), v_AddRobinEdgeContribution(), and v_AddRobinMassMatrix().

std::map<int, StdRegions::NormalVector> Nektar::LocalRegions::Expansion2D::m_edgeNormals
protected

Definition at line 123 of file Expansion2D.h.

Referenced by Nektar::LocalRegions::TriExp::v_ComputeEdgeNormal(), Nektar::LocalRegions::QuadExp::v_ComputeEdgeNormal(), Nektar::LocalRegions::NodalTriExp::v_ComputeEdgeNormal(), v_GetEdgeNormal(), and v_NegateEdgeNormal().

int Nektar::LocalRegions::Expansion2D::m_elementFaceLeft
protected

Definition at line 127 of file Expansion2D.h.

Referenced by Expansion2D(), GetLeftAdjacentElementFace(), and SetAdjacentElementExp().

int Nektar::LocalRegions::Expansion2D::m_elementFaceRight
protected

Definition at line 128 of file Expansion2D.h.

Referenced by Expansion2D(), GetRightAdjacentElementFace(), and SetAdjacentElementExp().

Expansion3DWeakPtr Nektar::LocalRegions::Expansion2D::m_elementLeft
protected

Definition at line 125 of file Expansion2D.h.

Referenced by GetLeftAdjacentElementExp(), GetRightAdjacentElementExp(), and SetAdjacentElementExp().

Expansion3DWeakPtr Nektar::LocalRegions::Expansion2D::m_elementRight
protected

Definition at line 126 of file Expansion2D.h.

Referenced by GetRightAdjacentElementExp(), and SetAdjacentElementExp().

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

Definition at line 124 of file Expansion2D.h.

Referenced by AddHDGHelmholtzEdgeTerms(), AddNormTraceInt(), v_AddEdgeNormBoundaryInt(), v_EdgeNormalNegated(), v_GenMatrix(), and v_NegateEdgeNormal().

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

Definition at line 122 of file Expansion2D.h.

Referenced by v_AddEdgeNormBoundaryInt().