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

#include <Expansion3D.h>

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

Public Member Functions

 Expansion3D (SpatialDomains::Geometry3DSharedPtr pGeom)
 
virtual ~Expansion3D ()
 
void SetTraceToGeomOrientation (Array< OneD, NekDouble > &inout)
 Align trace orientation with the geometry orientation. More...
 
void SetFaceToGeomOrientation (const int face, Array< OneD, NekDouble > &inout)
 Align face orientation with the geometry orientation. More...
 
void AddHDGHelmholtzFaceTerms (const NekDouble tau, const int edge, Array< OneD, NekDouble > &facePhys, const StdRegions::VarCoeffMap &dirForcing, Array< OneD, NekDouble > &outarray)
 
void AddNormTraceInt (const int dir, Array< OneD, ExpansionSharedPtr > &FaceExp, Array< OneD, Array< OneD, NekDouble > > &faceCoeffs, Array< OneD, NekDouble > &outarray)
 
void AddNormTraceInt (const int dir, Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &FaceExp, Array< OneD, NekDouble > &outarray, const StdRegions::VarCoeffMap &varcoeffs)
 
void AddFaceBoundaryInt (const int face, ExpansionSharedPtr &FaceExp, Array< OneD, NekDouble > &facePhys, Array< OneD, NekDouble > &outarray, const StdRegions::VarCoeffMap &varcoeffs=StdRegions::NullVarCoeffMap)
 
SpatialDomains::Geometry3DSharedPtr GetGeom3D () const
 
void v_ReOrientTracePhysMap (const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1)
 
void v_NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
Array< OneD, unsigned int > GetEdgeInverseBoundaryMap (int eid)
 
Array< OneD, unsigned int > GetTraceInverseBoundaryMap (int fid, StdRegions::Orientation faceOrient=StdRegions::eNoOrientation, int P1=-1, int P2=-1)
 
void GetInverseBoundaryMaps (Array< OneD, unsigned int > &vmap, Array< OneD, Array< OneD, unsigned int > > &emap, Array< OneD, Array< OneD, unsigned int > > &fmap)
 
- 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)
 
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)
 
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)
 
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) const
 
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)
 
virtual 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
 
- 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 (void) const
 
std::shared_ptr< StdExpansionGetLinStdExp (void) const
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion ()
 
bool IsNodalNonTensorialExp ()
 
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space. More...
 
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space. More...
 
void FwdTrans_BndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble Integral (const Array< OneD, const NekDouble > &inarray)
 This function integrates the specified function over the domain. More...
 
void FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 This function fills the array outarray with the mode-th mode of the expansion. More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 this function calculates the inner product of a given function f with the different modes of the expansion More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check)
 
void IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
int GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void SetElmtId (const int id)
 Set the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2=NullNekDouble1DArray, Array< OneD, NekDouble > &coords_3=NullNekDouble1DArray)
 this function returns the physical coordinates of the quadrature points of the expansion More...
 
void GetCoord (const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
 given the coordinates of a point of the element in the local collapsed coordinate system, this function calculates the physical coordinates of the point More...
 
DNekMatSharedPtr GetStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr GetStdStaticCondMatrix (const StdMatrixKey &mkey)
 
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 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, 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 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::StdExpansion3D
 StdExpansion3D ()
 
 StdExpansion3D (int numcoeffs, const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
 
 StdExpansion3D (const StdExpansion3D &T)
 
virtual ~StdExpansion3D ()
 
void PhysTensorDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d1, Array< OneD, NekDouble > &outarray_d2, Array< OneD, NekDouble > &outarray_d3)
 Calculate the 3D derivative in the local tensor/collapsed coordinate at the physical points. More...
 
void BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
 
void IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
 
int GetNedges () const
 return the number of edges in 3D expansion More...
 
int GetEdgeNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th edge. More...
 
void GetEdgeInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards)
 

Protected Member Functions

virtual void v_DGDeriv (const int dir, const Array< OneD, const NekDouble > &incoeffs, Array< OneD, ExpansionSharedPtr > &FaceExp, Array< OneD, Array< OneD, NekDouble > > &faceCoeffs, Array< OneD, NekDouble > &out_d)
 Evaluate coefficients of weak deriviative in the direction dir given the input coefficicents incoeffs and the imposed boundary values in EdgeExp (which will have its phys space updated). More...
 
virtual DNekMatSharedPtr v_GenMatrix (const StdRegions::StdMatrixKey &mkey)
 
virtual void v_AddFaceNormBoundaryInt (const int face, const ExpansionSharedPtr &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddRobinMassMatrix (const int face, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
virtual StdRegions::Orientation v_GetTraceOrient (int face)
 
virtual void v_GetTracePhysVals (const int face, const StdRegions::StdExpansionSharedPtr &FaceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient)
 Extract the physical values along face face from inarray into outarray following the local face orientation and point distribution defined by defined in FaceExp. More...
 
void GetPhysFaceVarCoeffsFromElement (const int face, ExpansionSharedPtr &FaceExp, const Array< OneD, const NekDouble > &varcoeff, Array< OneD, NekDouble > &outarray)
 
virtual Array< OneD, NekDoublev_GetnFacecdotMF (const int dir, const int face, ExpansionSharedPtr &FaceExp_f, const Array< OneD, const Array< OneD, NekDouble > > &normals, const StdRegions::VarCoeffMap &varcoeffs)
 
virtual DNekMatSharedPtr v_BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
 
virtual DNekMatSharedPtr v_BuildInverseTransformationMatrix (const DNekScalMatSharedPtr &transformationmatrix)
 Build inverse and inverse transposed transformation matrix: \(\mathbf{R^{-1}}\) and \(\mathbf{R^{-T}}\). More...
 
virtual DNekMatSharedPtr v_BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd)
 
virtual const NormalVectorv_GetTraceNormal (const int face) const
 
- Protected Member Functions inherited from Nektar::LocalRegions::Expansion
void ComputeLaplacianMetric ()
 
void ComputeQuadratureMetric ()
 
void ComputeGmatcdotMF (const Array< TwoD, const NekDouble > &df, const Array< OneD, const NekDouble > &direction, Array< OneD, Array< OneD, NekDouble > > &dfdir)
 
virtual void v_MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_DivideByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_ComputeLaplacianMetric ()
 
virtual void v_GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3)
 
Array< OneD, NekDoublev_GetMF (const int dir, const int shapedim, const StdRegions::VarCoeffMap &varcoeffs)
 
Array< OneD, NekDoublev_GetMFDiv (const int dir, const StdRegions::VarCoeffMap &varcoeffs)
 
Array< OneD, NekDoublev_GetMFMag (const int dir, const StdRegions::VarCoeffMap &varcoeffs)
 
virtual DNekScalMatSharedPtr v_GetLocMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual 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 NekDouble v_VectorFlux (const Array< OneD, Array< OneD, NekDouble > > &vec)
 
virtual void v_AlignVectorToCollapsedDir (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual void v_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_GetTraceQFactors (const int trace, Array< OneD, NekDouble > &outarray)
 
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 (void)
 
virtual void v_SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
virtual void v_SetUpPhysNormals (const int id)
 
virtual void v_AddRobinTraceContribution (const int traceid, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs)
 
- 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_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase_SumFac (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GeneralMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
 
void LaplacianMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp_MatFree (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void HelmholtzMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_SetCoeffsToOrientation (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>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals)
 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)
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion3D
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. More...
 
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 > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)=0
 
virtual void v_IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)=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)
 
virtual NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray)
 Integrates the specified function over the domain. More...
 
virtual int v_GetNedges (void) const
 
virtual int v_GetEdgeNcoeffs (const int i) const
 
virtual void v_GetEdgeInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards)
 
virtual void v_GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 

Protected Attributes

std::map< int, NormalVectorm_faceNormals
 
- Protected Attributes inherited from Nektar::LocalRegions::Expansion
LibUtilities::NekManager< IndexMapKey, IndexMapValues, IndexMapKey::opLessm_indexMapManager
 
std::vector< ExpansionWeakPtrm_traceExp
 
SpatialDomains::GeometrySharedPtr m_geom
 
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
 
MetricMap m_metrics
 
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 Attributes

std::vector< bool > m_requireNeg
 

Detailed Description

Definition at line 57 of file Expansion3D.h.

Constructor & Destructor Documentation

◆ Expansion3D()

Nektar::LocalRegions::Expansion3D::Expansion3D ( SpatialDomains::Geometry3DSharedPtr  pGeom)
inline

Definition at line 61 of file Expansion3D.h.

61 : Expansion(pGeom), StdExpansion3D(), m_requireNeg() {}
std::vector< bool > m_requireNeg
Definition: Expansion3D.h:188
Expansion(SpatialDomains::GeometrySharedPtr pGeom)
Definition: Expansion.cpp:47

◆ ~Expansion3D()

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

Definition at line 62 of file Expansion3D.h.

62 {}

Member Function Documentation

◆ AddFaceBoundaryInt()

void Nektar::LocalRegions::Expansion3D::AddFaceBoundaryInt ( const int  face,
ExpansionSharedPtr FaceExp,
Array< OneD, NekDouble > &  facePhys,
Array< OneD, NekDouble > &  outarray,
const StdRegions::VarCoeffMap varcoeffs = StdRegions::NullVarCoeffMap 
)
inline

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

Definition at line 342 of file Expansion3D.cpp.

348  {
349  boost::ignore_unused(varcoeffs);
350 
351  int i;
352  int order_f = FaceExp->GetNcoeffs();
353  Array<OneD, NekDouble> coeff(order_f);
354 
355  IndexMapKey ikey(eFaceToElement, DetShapeType(),
357  face, GetTraceOrient(face));
359 
360 // StdRegions::VarCoeffType VarCoeff[3] = {StdRegions::eVarCoeffD00,
361 // StdRegions::eVarCoeffD11,
362 // StdRegions::eVarCoeffD22};
363 // StdRegions::VarCoeffMap::const_iterator x;
364 // Array<OneD, NekDouble> varcoeff_work(nquad_e);
365 //
366 ///// @TODO Variable coeffs
367 // if ((x = varcoeffs.find(VarCoeff[0])) != varcoeffs.end())
368 // {
369 // GetPhysEdgeVarCoeffsFromElement(edge,EdgeExp,x->second,varcoeff_work);
370 // Vmath::Vmul(nquad_e,varcoeff_work,1,EdgeExp->GetPhys(),1,EdgeExp->UpdatePhys(),1);
371 // }
372 
373  FaceExp->IProductWRTBase(facePhys, coeff);
374 
375  // add data to out array
376  for(i = 0; i < order_f; ++i)
377  {
378  outarray[(*map)[i].index] += (*map)[i].sign*coeff[i];
379  }
380  }
StdRegions::Orientation GetTraceOrient(int trace)
Definition: Expansion.h:164
IndexMapValuesSharedPtr GetIndexMap(const IndexMapKey &ikey)
Definition: Expansion.h:138
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:376
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:171
std::shared_ptr< IndexMapValues > IndexMapValuesSharedPtr
Definition: IndexMapKey.h:126

References Nektar::LocalRegions::eFaceToElement.

◆ AddHDGHelmholtzFaceTerms()

void Nektar::LocalRegions::Expansion3D::AddHDGHelmholtzFaceTerms ( const NekDouble  tau,
const int  edge,
Array< OneD, NekDouble > &  facePhys,
const StdRegions::VarCoeffMap dirForcing,
Array< OneD, NekDouble > &  outarray 
)
inline

Definition at line 52 of file Expansion3D.cpp.

58  {
59  ExpansionSharedPtr FaceExp = GetTraceExp(face);
60  int i,j,n;
61  int nquad_f = FaceExp->GetNumPoints(0)*FaceExp->GetNumPoints(1);
62  int order_f = FaceExp->GetNcoeffs();
63  int coordim = GetCoordim();
64  int ncoeffs = GetNcoeffs();
65  bool mmf = (varcoeffs.find(StdRegions::eVarCoeffMF1x)
66  != varcoeffs.end());
67 
68  Array<OneD, NekDouble> inval (nquad_f);
69  Array<OneD, NekDouble> outcoeff(order_f);
70  Array<OneD, NekDouble> tmpcoeff(ncoeffs);
71 
72  const Array<OneD, const Array<OneD, NekDouble> > &normals
73  = GetTraceNormal(face);
74 
76 
77  DNekVec Coeffs(ncoeffs,outarray,eWrapper);
78  DNekVec Tmpcoeff(ncoeffs,tmpcoeff,eWrapper);
79 
80  IndexMapKey ikey(eFaceToElement, DetShapeType(),
82  face, GetTraceOrient(face));
84 
88 
89  // @TODO Variable coefficients
90  /*
91  StdRegions::VarCoeffType VarCoeff[3] = {StdRegions::eVarCoeffD00,
92  StdRegions::eVarCoeffD11,
93  StdRegions::eVarCoeffD22};
94  Array<OneD, NekDouble> varcoeff_work(nquad_f);
95  StdRegions::VarCoeffMap::const_iterator x;
96  ///// @TODO: What direction to use here??
97  if ((x = varcoeffs.find(VarCoeff[0])) != varcoeffs.end())
98  {
99  GetPhysFaceVarCoeffsFromElement(face,FaceExp,x->second,varcoeff_work);
100  Vmath::Vmul(nquad_f,varcoeff_work,1,FaceExp->GetPhys(),1,FaceExp->UpdatePhys(),1);
101  }
102  */
103 
104  //================================================================
105  // Add F = \tau <phi_i,in_phys>
106  // Fill face and take inner product
107  FaceExp->IProductWRTBase(facePhys, outcoeff);
108 
109  for(i = 0; i < order_f; ++i)
110  {
111  outarray[(*map)[i].index] += (*map)[i].sign*tau*outcoeff[i];
112  }
113  //================================================================
114 
115 
116  //===============================================================
117  // Add -\sum_i D_i^T M^{-1} G_i + E_i M^{-1} G_i =
118  // \sum_i D_i M^{-1} G_i term
119 
120  // Three independent direction
121  for(n = 0; n < coordim; ++n)
122  {
123  if (mmf) {
125  Weight[StdRegions::eVarCoeffMass] = v_GetMFMag(n,varcoeffs);
126 
127  MatrixKey invMasskey( StdRegions::eInvMass,
128  DetShapeType(), *this,
130  Weight);
131 
132  invMass = *GetLocMatrix(invMasskey);
133 
134  Array<OneD, NekDouble> ncdotMF_f =
135  v_GetnFacecdotMF(n, face, FaceExp, normals, varcoeffs);
136 
137  Vmath::Vmul(nquad_f, ncdotMF_f, 1, facePhys, 1, inval, 1);
138  }
139  else {
140  Vmath::Vmul(nquad_f, normals[n], 1, facePhys, 1, inval, 1);
141  }
142 
143  NekDouble scale = invMass.Scale();
144  const NekDouble *data = invMass.GetRawPtr();
145 
146  // @TODO Multiply by variable coefficients
147  // @TODO: Document this (probably not needed)
148  /*
149  StdRegions::VarCoeffMap::const_iterator x;
150  if ((x = varcoeffs.find(VarCoeff[n])) != varcoeffs.end())
151  {
152  GetPhysEdgeVarCoeffsFromElement(edge,FaceExp,x->second,varcoeff_work);
153  Vmath::Vmul(nquad_f,varcoeff_work,1,FaceExp->GetPhys(),1,FaceExp->UpdatePhys(),1);
154  }
155  */
156 
157  FaceExp->IProductWRTBase(inval, outcoeff);
158 
159  // M^{-1} G
160  for(i = 0; i < ncoeffs; ++i)
161  {
162  tmpcoeff[i] = 0;
163  for(j = 0; j < order_f; ++j)
164  {
165  tmpcoeff[i] += scale*data[i+(*map)[j].index*ncoeffs]*(*map)[j].sign*outcoeff[j];
166  }
167  }
168 
169  if (mmf)
170  {
171  StdRegions::VarCoeffMap VarCoeffDirDeriv;
172  VarCoeffDirDeriv[StdRegions::eVarCoeffMF]
173  = v_GetMF(n,coordim,varcoeffs);
174  VarCoeffDirDeriv[StdRegions::eVarCoeffMFDiv]
175  = v_GetMFDiv(n,varcoeffs);
176 
177  MatrixKey Dmatkey(StdRegions::eWeakDirectionalDeriv,
178  DetShapeType(), *this,
180  VarCoeffDirDeriv);
181 
182  DNekScalMat &Dmat = *GetLocMatrix(Dmatkey);
183 
184  Coeffs = Coeffs + Dmat*Tmpcoeff;
185  }
186  else
187  {
188  DNekScalMat &Dmat = *GetLocMatrix(DerivType[n]);
189  Coeffs = Coeffs + Dmat*Tmpcoeff;
190  }
191 
192  /*
193  if(varcoeffs.find(VarCoeff[n]) != varcoeffs.end())
194  {
195  MatrixKey mkey(DerivType[n], DetExpansionType(), *this, StdRegions::NullConstFactorMap, varcoeffs);
196  DNekScalMat &Dmat = *GetLocMatrix(mkey);
197  Coeffs = Coeffs + Dmat*Tmpcoeff;
198  }
199 
200  else
201  {
202  DNekScalMat &Dmat = *GetLocMatrix(DerivType[n]);
203  Coeffs = Coeffs + Dmat*Tmpcoeff;
204  }
205  */
206  }
207  }
virtual Array< OneD, NekDouble > v_GetnFacecdotMF(const int dir, const int face, ExpansionSharedPtr &FaceExp_f, const Array< OneD, const Array< OneD, NekDouble > > &normals, const StdRegions::VarCoeffMap &varcoeffs)
Array< OneD, NekDouble > v_GetMF(const int dir, const int shapedim, const StdRegions::VarCoeffMap &varcoeffs)
Definition: Expansion.cpp:434
const NormalVector & GetTraceNormal(const int id) const
Definition: Expansion.h:222
Array< OneD, NekDouble > v_GetMFMag(const int dir, const StdRegions::VarCoeffMap &varcoeffs)
Definition: Expansion.cpp:531
Array< OneD, NekDouble > v_GetMFDiv(const int dir, const StdRegions::VarCoeffMap &varcoeffs)
Definition: Expansion.cpp:498
ExpansionSharedPtr GetTraceExp(const int traceid)
Definition: Expansion.h:417
DNekScalMatSharedPtr GetLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:90
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:124
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:272
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:315
NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag > DNekScalMat
NekVector< NekDouble > DNekVec
Definition: NekTypeDefs.hpp:48
double NekDouble
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:192

References Nektar::LocalRegions::eFaceToElement, Nektar::StdRegions::eInvMass, 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::NullConstFactorMap, and Vmath::Vmul().

◆ AddNormTraceInt() [1/2]

void Nektar::LocalRegions::Expansion3D::AddNormTraceInt ( const int  dir,
Array< OneD, const NekDouble > &  inarray,
Array< OneD, ExpansionSharedPtr > &  FaceExp,
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 242 of file Expansion3D.cpp.

248  {
249  int i,f,cnt;
250  int order_f,nquad_f;
251  int nfaces = GetNtraces();
252 
253  cnt = 0;
254  for(f = 0; f < nfaces; ++f)
255  {
256  order_f = FaceExp[f]->GetNcoeffs();
257  nquad_f = FaceExp[f]->GetNumPoints(0)*FaceExp[f]->GetNumPoints(1);
258 
259  const Array<OneD, const Array<OneD, NekDouble> > &normals = GetTraceNormal(f);
260  Array<OneD, NekDouble> faceCoeffs(order_f);
261  Array<OneD, NekDouble> facePhys (nquad_f);
262 
263  for(i = 0; i < order_f; ++i)
264  {
265  faceCoeffs[i] = inarray[i+cnt];
266  }
267  cnt += order_f;
268 
269  FaceExp[f]->BwdTrans(faceCoeffs, facePhys);
270 
271  // Multiply by variable coefficient
272  /// @TODO: Document this
273 // StdRegions::VarCoeffType VarCoeff[3] = {StdRegions::eVarCoeffD00,
274 // StdRegions::eVarCoeffD11,
275 // StdRegions::eVarCoeffD22};
276 // StdRegions::VarCoeffMap::const_iterator x;
277 // Array<OneD, NekDouble> varcoeff_work(nquad_e);
278 
279 // if ((x = varcoeffs.find(VarCoeff[dir])) != varcoeffs.end())
280 // {
281 // GetPhysEdgeVarCoeffsFromElement(e,EdgeExp[e],x->second,varcoeff_work);
282 // Vmath::Vmul(nquad_e,varcoeff_work,1,EdgeExp[e]->GetPhys(),1,EdgeExp[e]->UpdatePhys(),1);
283 // }
284  StdRegions::VarCoeffMap::const_iterator x;
285  if ((x = varcoeffs.find(StdRegions::eVarCoeffMF1x))
286  != varcoeffs.end())
287  {
288  Array<OneD, NekDouble> ncdotMF_f =
289  v_GetnFacecdotMF(dir, f, FaceExp[f], normals,
290  varcoeffs);
291 
292  Vmath::Vmul(nquad_f, ncdotMF_f, 1,
293  facePhys, 1,
294  facePhys, 1);
295  }
296  else
297  {
298  Vmath::Vmul(nquad_f, normals[dir], 1,
299  facePhys, 1,
300  facePhys, 1);
301  }
302 
303  AddFaceBoundaryInt(f, FaceExp[f], facePhys, outarray,
304  varcoeffs);
305  }
306  }
void AddFaceBoundaryInt(const int face, ExpansionSharedPtr &FaceExp, Array< OneD, NekDouble > &facePhys, Array< OneD, NekDouble > &outarray, const StdRegions::VarCoeffMap &varcoeffs=StdRegions::NullVarCoeffMap)
int GetNtraces() const
Returns the number of trace elements connected to this element.
Definition: StdExpansion.h:360

References Nektar::StdRegions::eVarCoeffMF1x, and Vmath::Vmul().

◆ AddNormTraceInt() [2/2]

void Nektar::LocalRegions::Expansion3D::AddNormTraceInt ( const int  dir,
Array< OneD, ExpansionSharedPtr > &  FaceExp,
Array< OneD, Array< OneD, NekDouble > > &  faceCoeffs,
Array< OneD, NekDouble > &  outarray 
)
inline

Definition at line 309 of file Expansion3D.cpp.

314  {
315  int f, cnt;
316  int order_f, nquad_f;
317  int nfaces = GetNtraces();
318 
319  cnt = 0;
320  for(f = 0; f < nfaces; ++f)
321  {
322  order_f = FaceExp[f]->GetNcoeffs();
323  nquad_f = FaceExp[f]->GetNumPoints(0)*FaceExp[f]->GetNumPoints(1);
324 
325  const Array<OneD, const Array<OneD, NekDouble> > &normals =
326  GetTraceNormal(f);
327  Array<OneD, NekDouble> facePhys(nquad_f);
328 
329  cnt += order_f;
330 
331  FaceExp[f]->BwdTrans(faceCoeffs[f], facePhys);
332 
333  Vmath::Vmul(nquad_f, normals[dir], 1, facePhys, 1, facePhys, 1);
334 
335  AddFaceBoundaryInt(f, FaceExp[f], facePhys, outarray);
336  }
337  }

References Vmath::Vmul().

◆ GetEdgeInverseBoundaryMap()

Array< OneD, unsigned int > Nektar::LocalRegions::Expansion3D::GetEdgeInverseBoundaryMap ( int  eid)

Definition at line 1990 of file Expansion3D.cpp.

1992  {
1993  int n, j;
1994  int nEdgeCoeffs;
1995  int nBndCoeffs = NumBndryCoeffs();
1996 
1997  Array<OneD, unsigned int> bmap(nBndCoeffs);
1998  GetBoundaryMap(bmap);
1999 
2000  // Map from full system to statically condensed system (i.e reverse
2001  // GetBoundaryMap)
2002  map<int, int> invmap;
2003  for (j = 0; j < nBndCoeffs; ++j)
2004  {
2005  invmap[bmap[j]] = j;
2006  }
2007 
2008  // Number of interior edge coefficients
2009  nEdgeCoeffs = GetEdgeNcoeffs(eid) - 2;
2010 
2012 
2013  Array<OneD, unsigned int> edgemaparray(nEdgeCoeffs);
2014  StdRegions::Orientation eOrient =
2015  geom->GetEorient(eid);
2016  Array<OneD, unsigned int> maparray =
2017  Array<OneD, unsigned int>(nEdgeCoeffs);
2018  Array<OneD, int> signarray =
2019  Array<OneD, int>(nEdgeCoeffs, 1);
2020 
2021  // maparray is the location of the edge within the matrix
2022  GetEdgeInteriorToElementMap(eid, maparray, signarray, eOrient);
2023 
2024  for (n = 0; n < nEdgeCoeffs; ++n)
2025  {
2026  edgemaparray[n] = invmap[maparray[n]];
2027  }
2028 
2029  return edgemaparray;
2030  }
SpatialDomains::Geometry3DSharedPtr GetGeom3D() const
Definition: Expansion3D.h:191
int GetEdgeNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th edge.
void GetEdgeInteriorToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards)
void GetBoundaryMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:687
std::shared_ptr< Geometry3D > Geometry3DSharedPtr
Definition: Geometry3D.h:52

◆ GetGeom3D()

SpatialDomains::Geometry3DSharedPtr Nektar::LocalRegions::Expansion3D::GetGeom3D ( ) const
inline

Definition at line 191 of file Expansion3D.h.

192  {
193  return std::dynamic_pointer_cast<SpatialDomains::Geometry3D>(m_geom);
194  }
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:272

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

Referenced by Nektar::LocalRegions::HexExp::v_GetFaceDGForwards().

◆ GetInverseBoundaryMaps()

void Nektar::LocalRegions::Expansion3D::GetInverseBoundaryMaps ( Array< OneD, unsigned int > &  vmap,
Array< OneD, Array< OneD, unsigned int > > &  emap,
Array< OneD, Array< OneD, unsigned int > > &  fmap 
)

Definition at line 2148 of file Expansion3D.cpp.

2152  {
2153  int n, j;
2154  int nEdgeCoeffs;
2155  int nFaceCoeffs;
2156 
2157  int nBndCoeffs = NumBndryCoeffs();
2158 
2159  Array<OneD, unsigned int> bmap(nBndCoeffs);
2160  GetBoundaryMap(bmap);
2161 
2162  // Map from full system to statically condensed system (i.e reverse
2163  // GetBoundaryMap)
2164  map<int, int> reversemap;
2165  for (j = 0; j < bmap.size(); ++j)
2166  {
2167  reversemap[bmap[j]] = j;
2168  }
2169 
2170  int nverts = GetNverts();
2171  vmap = Array<OneD, unsigned int>(nverts);
2172  for (n = 0; n < nverts; ++n)
2173  {
2174  int id = GetVertexMap(n);
2175  vmap[n] = reversemap[id]; // not sure what should be true here.
2176  }
2177 
2178  int nedges = GetNedges();
2179  emap = Array<OneD, Array<OneD, unsigned int> >(nedges);
2180 
2181  for(int eid = 0; eid < nedges; ++eid)
2182  {
2183  // Number of interior edge coefficients
2184  nEdgeCoeffs = GetEdgeNcoeffs(eid) - 2;
2185 
2186  Array<OneD, unsigned int> edgemaparray(nEdgeCoeffs);
2187  Array<OneD, unsigned int> maparray =
2188  Array<OneD, unsigned int>(nEdgeCoeffs);
2189  Array<OneD, int> signarray =
2190  Array<OneD, int>(nEdgeCoeffs, 1);
2191 
2192  // maparray is the location of the edge within the matrix
2193  GetEdgeInteriorToElementMap(eid, maparray, signarray, StdRegions::eForwards);
2194 
2195  for (n = 0; n < nEdgeCoeffs; ++n)
2196  {
2197  edgemaparray[n] = reversemap[maparray[n]];
2198  }
2199  emap[eid] = edgemaparray;
2200  }
2201 
2202  int nfaces = GetNtraces();
2203  fmap = Array<OneD, Array<OneD, unsigned int> >(nfaces);
2204 
2205  for(int fid = 0; fid < nfaces; ++fid)
2206  {
2207  // Number of interior face coefficients
2208  nFaceCoeffs = GetTraceIntNcoeffs(fid);
2209 
2210  Array<OneD, unsigned int> facemaparray(nFaceCoeffs);
2211  Array<OneD, unsigned int> maparray =
2212  Array<OneD, unsigned int>(nFaceCoeffs);
2213  Array<OneD, int> signarray =
2214  Array<OneD, int>(nFaceCoeffs, 1);
2215 
2216  // maparray is the location of the face within the matrix
2217  GetTraceInteriorToElementMap(fid,maparray, signarray,
2219 
2220  for (n = 0; n < nFaceCoeffs; ++n)
2221  {
2222  facemaparray[n] = reversemap[maparray[n]];
2223  }
2224 
2225  fmap[fid] = facemaparray;
2226  }
2227  }
int GetNedges() const
return the number of edges in 3D expansion
int GetTraceIntNcoeffs(const int i) const
Definition: StdExpansion.h:270
int GetVertexMap(const int localVertexId, bool useCoeffPacking=false)
Definition: StdExpansion.h:697
int GetNverts() const
This function returns the number of vertices of the expansion domain.
Definition: StdExpansion.h:249
void GetTraceInteriorToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eForwards)
Definition: StdExpansion.h:714

References Nektar::StdRegions::eDir1FwdDir1_Dir2FwdDir2, and Nektar::StdRegions::eForwards.

◆ GetPhysFaceVarCoeffsFromElement()

void Nektar::LocalRegions::Expansion3D::GetPhysFaceVarCoeffsFromElement ( const int  face,
ExpansionSharedPtr FaceExp,
const Array< OneD, const NekDouble > &  varcoeff,
Array< OneD, NekDouble > &  outarray 
)
protected

Definition at line 210 of file Expansion3D.cpp.

215  {
216  Array<OneD, NekDouble> tmp(GetNcoeffs());
217  Array<OneD, NekDouble> facetmp(FaceExp->GetNcoeffs());
218 
219  // FwdTrans varcoeffs
220  FwdTrans(varcoeff, tmp);
221 
222  // Map to edge
223  Array<OneD, unsigned int> emap;
224  Array<OneD, int> sign;
225 
226  GetTraceToElementMap(face, emap, sign, GetTraceOrient(face));
227 
228  for (unsigned int i = 0; i < FaceExp->GetNcoeffs(); ++i)
229  {
230  facetmp[i] = tmp[emap[i]];
231  }
232 
233  // BwdTrans
234  FaceExp->BwdTrans(facetmp, outarray);
235  }
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:15
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:703
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 sign.

◆ GetTraceInverseBoundaryMap()

Array< OneD, unsigned int > Nektar::LocalRegions::Expansion3D::GetTraceInverseBoundaryMap ( int  fid,
StdRegions::Orientation  faceOrient = StdRegions::eNoOrientation,
int  P1 = -1,
int  P2 = -1 
)

Definition at line 2032 of file Expansion3D.cpp.

2037  {
2038  int n,j;
2039  int nFaceCoeffs;
2040 
2041  int nBndCoeffs = NumBndryCoeffs();
2042 
2043  Array<OneD, unsigned int> bmap(nBndCoeffs);
2044  GetBoundaryMap(bmap);
2045 
2046  // Map from full system to statically condensed system (i.e reverse
2047  // GetBoundaryMap)
2048  map<int, int> reversemap;
2049  for (j = 0; j < bmap.size(); ++j)
2050  {
2051  reversemap[bmap[j]] = j;
2052  }
2053 
2054  // Number of interior face coefficients
2055  nFaceCoeffs = GetTraceIntNcoeffs(fid);
2056 
2057  StdRegions::Orientation fOrient;
2058  Array<OneD, unsigned int> maparray =
2059  Array<OneD, unsigned int>(nFaceCoeffs);
2060  Array<OneD, int> signarray =
2061  Array<OneD, int>(nFaceCoeffs, 1);
2062 
2063  if(faceOrient == StdRegions::eNoOrientation)
2064  {
2065  fOrient = GetTraceOrient(fid);
2066  }
2067  else
2068  {
2069  fOrient = faceOrient;
2070  }
2071 
2072  // maparray is the location of the face within the matrix
2073  GetTraceInteriorToElementMap(fid, maparray, signarray, fOrient);
2074 
2075  Array<OneD, unsigned int> facemaparray;
2076  int locP1,locP2;
2077  GetTraceNumModes(fid,locP1,locP2,fOrient);
2078 
2079  if(P1 == -1)
2080  {
2081  P1 = locP1;
2082  }
2083  else
2084  {
2085  ASSERTL1(P1 <= locP1,"Expect value of passed P1 to "
2086  "be lower or equal to face num modes");
2087  }
2088 
2089  if(P2 == -1)
2090  {
2091  P2 = locP2;
2092  }
2093  else
2094  {
2095  ASSERTL1(P2 <= locP2,"Expect value of passed P2 to "
2096  "be lower or equal to face num modes");
2097  }
2098 
2099  switch(GetGeom3D()->GetFace(fid)->GetShapeType())
2100  {
2102  {
2103  if(((P1-3)>0)&&((P2-3)>0))
2104  {
2105  facemaparray= Array<OneD, unsigned int>(LibUtilities::StdTriData::getNumberOfCoefficients(P1-3,P2-3));
2106  int cnt = 0;
2107  int cnt1 = 0;
2108  for(n = 0; n < P1-3; ++n)
2109  {
2110  for(int m = 0; m < P2-3-n; ++m, ++cnt)
2111  {
2112  facemaparray[cnt] = reversemap[maparray[cnt1+m]];
2113  }
2114  cnt1 += locP2-3-n;
2115  }
2116  }
2117  }
2118  break;
2120  {
2121  if(((P1-2)>0)&&((P2-2)>0))
2122  {
2123  facemaparray = Array<OneD, unsigned int>(LibUtilities::StdQuadData::getNumberOfCoefficients(P1-2,P2-2));
2124  int cnt = 0;
2125  int cnt1 = 0;
2126  for(n = 0; n < P2-2; ++n)
2127  {
2128  for(int m = 0; m < P1-2; ++m, ++cnt)
2129  {
2130  facemaparray[cnt] = reversemap[maparray[cnt1+m]];
2131  }
2132  cnt1 += locP1-2;
2133  }
2134  }
2135  }
2136  break;
2137  default:
2138  {
2139  ASSERTL0(false, "Invalid shape type.");
2140  }
2141  break;
2142  }
2143 
2144 
2145  return facemaparray;
2146  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
void GetTraceNumModes(const int tid, int &numModes0, int &numModes1, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2)
Definition: StdExpansion.h:724
int getNumberOfCoefficients(int Na, int Nb)
Definition: ShapeType.hpp:137
int getNumberOfCoefficients(int Na, int Nb)
Definition: ShapeType.hpp:113

References ASSERTL0, ASSERTL1, Nektar::StdRegions::eNoOrientation, Nektar::LibUtilities::eQuadrilateral, Nektar::LibUtilities::eTriangle, and Nektar::LibUtilities::StdSegData::getNumberOfCoefficients().

◆ SetFaceToGeomOrientation()

void Nektar::LocalRegions::Expansion3D::SetFaceToGeomOrientation ( const int  face,
Array< OneD, NekDouble > &  inout 
)

Align face orientation with the geometry orientation.

Definition at line 385 of file Expansion3D.cpp.

387  {
388  int j,k;
389  int nface = GetTraceNcoeffs(face);
390  Array<OneD, NekDouble> f_in(nface);
391  Vmath::Vcopy(nface,&inout[0],1,&f_in[0],1);
392 
393  // retreiving face to element map for standard face orientation and
394  // for actual face orientation
395  IndexMapKey ikey1(eFaceToElement, DetShapeType(),
398  IndexMapValuesSharedPtr map1 = GetIndexMap(ikey1);
399  IndexMapKey ikey2(eFaceToElement, DetShapeType(),
401  face, GetTraceOrient(face));
402  IndexMapValuesSharedPtr map2 = GetIndexMap(ikey2);
403 
404  ASSERTL1((*map1).size() == (*map2).size(),
405  "There is an error with the GetTraceToElementMap");
406 
407  for(j = 0; j < (*map1).size(); ++j)
408  {
409  // j = index in the standard orientation
410  for(k = 0; k < (*map2).size(); ++k)
411  {
412  // k = index in the actual orientation
413  if((*map1)[j].index == (*map2)[k].index && k != j)
414  {
415  inout[k] = f_in[j];
416  //checking if sign is changing
417  if((*map1)[j].sign != (*map2)[k].sign)
418  inout[k] *= -1.0;
419  break;
420  }
421  }
422  }
423  }
int GetTraceNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th trace.
Definition: StdExpansion.h:265
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199

References ASSERTL1, Nektar::StdRegions::eDir1FwdDir1_Dir2FwdDir2, Nektar::LocalRegions::eFaceToElement, sign, and Vmath::Vcopy().

◆ SetTraceToGeomOrientation()

void Nektar::LocalRegions::Expansion3D::SetTraceToGeomOrientation ( Array< OneD, NekDouble > &  inout)

Align trace orientation with the geometry orientation.

Definition at line 428 of file Expansion3D.cpp.

429  {
430  int i,cnt = 0;
431  int nfaces = GetNtraces();
432 
433  Array<OneD, NekDouble> f_tmp;
434 
435  for(i = 0; i < nfaces; ++i)
436  {
437  SetFaceToGeomOrientation(i, f_tmp = inout + cnt);
438  cnt += GetTraceNcoeffs(i);
439  }
440  }
void SetFaceToGeomOrientation(const int face, Array< OneD, NekDouble > &inout)
Align face orientation with the geometry orientation.

◆ v_AddFaceNormBoundaryInt()

void Nektar::LocalRegions::Expansion3D::v_AddFaceNormBoundaryInt ( const int  face,
const ExpansionSharedPtr FaceExp,
const Array< OneD, const NekDouble > &  Fn,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Definition at line 1092 of file Expansion3D.cpp.

1097  {
1098  int i, j;
1099 
1100  /*
1101  * Coming into this routine, the velocity V will have been
1102  * multiplied by the trace normals to give the input vector Vn. By
1103  * convention, these normals are inwards facing for elements which
1104  * have FaceExp as their right-adjacent face. This conditional
1105  * statement therefore determines whether the normals must be
1106  * negated, since the integral being performed here requires an
1107  * outwards facing normal.
1108  */
1109  if (m_requireNeg.size() == 0)
1110  {
1111  m_requireNeg.resize(GetNtraces());
1112 
1113  for (i = 0; i < GetNtraces(); ++i)
1114  {
1115  m_requireNeg[i] = false;
1116 
1117  ExpansionSharedPtr faceExp = m_traceExp[i].lock();
1118 
1119  if (faceExp->GetRightAdjacentElementExp())
1120  {
1121  if (faceExp->GetRightAdjacentElementExp()->GetGeom()
1122  ->GetGlobalID() == GetGeom()->GetGlobalID())
1123  {
1124  m_requireNeg[i] = true;
1125  }
1126  }
1127  }
1128  }
1129 
1130  IndexMapKey ikey(eFaceToElement, DetShapeType(),
1132  face, GetTraceOrient(face));
1134 
1135  int order_e = (*map).size(); // Order of the element
1136  int n_coeffs = FaceExp->GetNcoeffs();
1137 
1138  Array<OneD, NekDouble> faceCoeffs(n_coeffs);
1139 
1140  if (n_coeffs != order_e) // Going to orthogonal space
1141  {
1142  Array<OneD, NekDouble> coeff(n_coeffs);
1143  Array<OneD, NekDouble> array(n_coeffs);
1144 
1145  FaceExp->FwdTrans(Fn, faceCoeffs);
1146 
1147  int NumModesElementMax = FaceExp->GetBasis(0)->GetNumModes();
1148  int NumModesElementMin = m_base[0]->GetNumModes();
1149 
1150  FaceExp->ReduceOrderCoeffs(NumModesElementMin,
1151  faceCoeffs,
1152  faceCoeffs);
1153 
1154  StdRegions::StdMatrixKey masskey(StdRegions::eMass,
1155  FaceExp->DetShapeType(),
1156  *FaceExp);
1157  FaceExp->MassMatrixOp(faceCoeffs, faceCoeffs, masskey);
1158 
1159  // Reorder coefficients for the lower degree face.
1160  int offset1 = 0, offset2 = 0;
1161 
1162  if (FaceExp->DetShapeType() == LibUtilities::eQuadrilateral)
1163  {
1164  for (i = 0; i < NumModesElementMin; ++i)
1165  {
1166  for (j = 0; j < NumModesElementMin; ++j)
1167  {
1168  faceCoeffs[offset1+j] =
1169  faceCoeffs[offset2+j];
1170  }
1171  offset1 += NumModesElementMin;
1172  offset2 += NumModesElementMax;
1173  }
1174 
1175  // Extract lower degree modes. TODO: Check this is correct.
1176  for (i = NumModesElementMin; i < NumModesElementMax; ++i)
1177  {
1178  for (j = NumModesElementMin; j < NumModesElementMax; ++j)
1179  {
1180  faceCoeffs[i*NumModesElementMax+j] = 0.0;
1181  }
1182  }
1183  }
1184 
1185  if (FaceExp->DetShapeType() == LibUtilities::eTriangle)
1186  {
1187 
1188  // Reorder coefficients for the lower degree face.
1189  int offset1 = 0, offset2 = 0;
1190 
1191  for (i = 0; i < NumModesElementMin; ++i)
1192  {
1193  for (j = 0; j < NumModesElementMin-i; ++j)
1194  {
1195  faceCoeffs[offset1+j] =
1196  faceCoeffs[offset2+j];
1197  }
1198  offset1 += NumModesElementMin-i;
1199  offset2 += NumModesElementMax-i;
1200  }
1201  }
1202 
1203  }
1204  else
1205  {
1206  FaceExp->IProductWRTBase(Fn, faceCoeffs);
1207  }
1208 
1209  if (m_requireNeg[face])
1210  {
1211  for (i = 0; i < order_e; ++i)
1212  {
1213  outarray[(*map)[i].index] -= (*map)[i].sign * faceCoeffs[i];
1214  }
1215  }
1216  else
1217  {
1218  for (i = 0; i < order_e; ++i)
1219  {
1220  outarray[(*map)[i].index] += (*map)[i].sign * faceCoeffs[i];
1221  }
1222  }
1223  }
SpatialDomains::GeometrySharedPtr GetGeom() const
Definition: Expansion.cpp:172
std::vector< ExpansionWeakPtr > m_traceExp
Definition: Expansion.h:271
Array< OneD, LibUtilities::BasisSharedPtr > m_base

References Nektar::LocalRegions::eFaceToElement, Nektar::StdRegions::eMass, Nektar::LibUtilities::eQuadrilateral, and Nektar::LibUtilities::eTriangle.

◆ v_AddRobinMassMatrix()

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

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1261 of file Expansion3D.cpp.

1265  {
1267  "Not set up for non boundary-interior expansions");
1268  ASSERTL1(inoutmat->GetRows() == inoutmat->GetColumns(),
1269  "Assuming that input matrix was square");
1270 
1271  int i,j;
1272  int id1,id2;
1273  ExpansionSharedPtr faceExp = m_traceExp[face].lock();
1274  int order_f = faceExp->GetNcoeffs();
1275 
1276  Array<OneD, unsigned int> map;
1277  Array<OneD, int> sign;
1278 
1279  StdRegions::VarCoeffMap varcoeffs;
1280  varcoeffs[StdRegions::eVarCoeffMass] = primCoeffs;
1281 
1282  LibUtilities::ShapeType shapeType =
1283  faceExp->DetShapeType();
1284 
1285  LocalRegions::MatrixKey mkey(StdRegions::eMass,
1286  shapeType,
1287  *faceExp,
1289  varcoeffs);
1290 
1291  DNekScalMat &facemat = *faceExp->GetLocMatrix(mkey);
1292 
1293  // Now need to identify a map which takes the local face
1294  // mass matrix to the matrix stored in inoutmat;
1295  // This can currently be deduced from the size of the matrix
1296 
1297  // - if inoutmat.m_rows() == v_NCoeffs() it is a full
1298  // matrix system
1299 
1300  // - if inoutmat.m_rows() == v_GetNverts() it is a vertex space
1301  // preconditioner.
1302 
1303  // - if inoutmat.m_rows() == v_NumBndCoeffs() it is a
1304  // boundary CG system
1305 
1306  // - if inoutmat.m_rows() == v_NumDGBndCoeffs() it is a
1307  // trace DG system; still needs implementing.
1308  int rows = inoutmat->GetRows();
1309 
1310  if (rows == GetNcoeffs())
1311  {
1312  GetTraceToElementMap(face,map,sign,GetTraceOrient(face));
1313  }
1314  else if (rows == GetNverts())
1315  {
1316  int nfvert = faceExp->GetNverts();
1317 
1318  // Need to find where linear vertices are in facemat
1319  Array<OneD, unsigned int> linmap;
1320  Array<OneD, int> linsign;
1321 
1322  // Use a linear expansion to get correct mapping
1323  GetLinStdExp()->GetTraceToElementMap(face,linmap,linsign,
1324  GetTraceOrient(face));
1325 
1326  // zero out sign map to remove all other modes
1327  sign = Array<OneD, int> (order_f,0);
1328  map = Array<OneD, unsigned int>(order_f,(unsigned int)0);
1329 
1330  int fmap;
1331  // Reset sign map to only have contribution from vertices
1332  for(i = 0; i < nfvert; ++i)
1333  {
1334  fmap = faceExp->GetVertexMap(i,true);
1335  sign[fmap] = 1;
1336 
1337  // need to reset map
1338  map[fmap] = linmap[i];
1339  }
1340  }
1341  else if(rows == NumBndryCoeffs())
1342  {
1343  int nbndry = NumBndryCoeffs();
1344  Array<OneD,unsigned int> bmap(nbndry);
1345 
1346  GetTraceToElementMap(face,map,sign,GetTraceOrient(face));
1347  GetBoundaryMap(bmap);
1348 
1349  for(i = 0; i < order_f; ++i)
1350  {
1351  for(j = 0; j < nbndry; ++j)
1352  {
1353  if(map[i] == bmap[j])
1354  {
1355  map[i] = j;
1356  break;
1357  }
1358  }
1359  ASSERTL1(j != nbndry,"Did not find number in map");
1360  }
1361  }
1362  else if (rows == NumDGBndryCoeffs())
1363  {
1364  // possibly this should be a separate method
1365  int cnt = 0;
1366  map = Array<OneD, unsigned int> (order_f);
1367  sign = Array<OneD, int> (order_f,1);
1368 
1369  IndexMapKey ikey1(eFaceToElement, DetShapeType(),
1371  GetBasisNumModes(2), face, GetTraceOrient(face));
1372  IndexMapValuesSharedPtr map1 = GetIndexMap(ikey1);
1373  IndexMapKey ikey2(eFaceToElement, DetShapeType(),
1375  GetBasisNumModes(2), face,
1377  IndexMapValuesSharedPtr map2 = GetIndexMap(ikey2);
1378 
1379  ASSERTL1((*map1).size() == (*map2).size(),
1380  "There is an error with the GetTraceToElementMap");
1381 
1382  for (i = 0; i < face; ++i)
1383  {
1384  cnt += GetTraceNcoeffs(i);
1385  }
1386 
1387  for(i = 0; i < (*map1).size(); ++i)
1388  {
1389  int idx = -1;
1390 
1391  for(j = 0; j < (*map2).size(); ++j)
1392  {
1393  if((*map1)[i].index == (*map2)[j].index)
1394  {
1395  idx = j;
1396  break;
1397  }
1398  }
1399 
1400  ASSERTL2(idx >= 0, "Index not found");
1401  map [i] = idx + cnt;
1402  sign[i] = (*map2)[idx].sign;
1403  }
1404  }
1405  else
1406  {
1407  ASSERTL0(false,"Could not identify matrix type from dimension");
1408  }
1409 
1410  for(i = 0; i < order_f; ++i)
1411  {
1412  id1 = map[i];
1413  for(j = 0; j < order_f; ++j)
1414  {
1415  id2 = map[j];
1416  (*inoutmat)(id1,id2) += facemat(i,j)*sign[i]*sign[j];
1417  }
1418  }
1419  }
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:274
std::shared_ptr< StdExpansion > GetLinStdExp(void) const
Definition: StdExpansion.h:386

References ASSERTL0, ASSERTL1, ASSERTL2, Nektar::StdRegions::eDir1FwdDir1_Dir2FwdDir2, Nektar::LocalRegions::eFaceToElement, Nektar::StdRegions::eMass, Nektar::StdRegions::eVarCoeffMass, Nektar::StdRegions::NullConstFactorMap, and sign.

◆ v_BuildInverseTransformationMatrix()

DNekMatSharedPtr Nektar::LocalRegions::Expansion3D::v_BuildInverseTransformationMatrix ( const DNekScalMatSharedPtr transformationmatrix)
protectedvirtual

Build inverse and inverse transposed transformation matrix: \(\mathbf{R^{-1}}\) and \(\mathbf{R^{-T}}\).

\f\mathbf{R^{-T}}=[\left[\begin{array}{ccc} \mathbf{I} & -\mathbf{R}_{ef} & -\mathbf{R}_{ve}+\mathbf{R}_{ve}\mathbf{R}_{vf} \ 0 & \mathbf{I} & \mathbf{R}_{ef} \ 0 & 0 & \mathbf{I}} \end{array}\right]\f]

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1863 of file Expansion3D.cpp.

1865  {
1866  int i, j, n, eid = 0, fid = 0;
1867  int nCoeffs = NumBndryCoeffs();
1868  NekDouble MatrixValue;
1869  DNekScalMat &R = (*transformationmatrix);
1870 
1871  // Define storage for vertex transpose matrix and zero all entries
1872  MatrixStorage storage = eFULL;
1873 
1874  // Inverse transformation matrix
1875  DNekMatSharedPtr inversetransformationmatrix =
1877  nCoeffs, nCoeffs, 0.0, storage);
1878  DNekMat &InvR = (*inversetransformationmatrix);
1879 
1880  int nVerts = GetNverts();
1881  int nEdges = GetNedges();
1882  int nFaces = GetNtraces();
1883 
1884  int nedgemodes = 0;
1885  int nfacemodes = 0;
1886  int nedgemodestotal = 0;
1887  int nfacemodestotal = 0;
1888 
1889  for (eid = 0; eid < nEdges; ++eid)
1890  {
1891  nedgemodes = GetEdgeNcoeffs(eid) - 2;
1892  nedgemodestotal += nedgemodes;
1893  }
1894 
1895  for (fid = 0; fid < nFaces; ++fid)
1896  {
1897  nfacemodes = GetTraceIntNcoeffs(fid);
1898  nfacemodestotal += nfacemodes;
1899  }
1900 
1901  Array<OneD, unsigned int>
1902  edgemodearray(nedgemodestotal);
1903  Array<OneD, unsigned int>
1904  facemodearray(nfacemodestotal);
1905 
1906  int offset = 0;
1907 
1908  // Create array of edge modes
1909  for (eid = 0; eid < nEdges; ++eid)
1910  {
1911  Array<OneD, unsigned int> edgearray
1913  nedgemodes = GetEdgeNcoeffs(eid) - 2;
1914 
1915  // Only copy if there are edge modes
1916  if (nedgemodes)
1917  {
1918  Vmath::Vcopy(nedgemodes, &edgearray[0], 1,
1919  &edgemodearray[offset], 1);
1920  }
1921 
1922  offset += nedgemodes;
1923  }
1924 
1925  offset = 0;
1926 
1927  // Create array of face modes
1928  for (fid = 0; fid < nFaces; ++fid)
1929  {
1930  Array<OneD, unsigned int> facearray
1932  nfacemodes = GetTraceIntNcoeffs(fid);
1933 
1934  // Only copy if there are face modes
1935  if (nfacemodes)
1936  {
1937  Vmath::Vcopy(nfacemodes, &facearray[0], 1,
1938  &facemodearray[offset], 1);
1939  }
1940 
1941  offset += nfacemodes;
1942  }
1943 
1944  // Vertex-edge/face
1945  for (i = 0; i < nVerts; ++i)
1946  {
1947  for (j = 0; j < nedgemodestotal; ++j)
1948  {
1949  InvR.SetValue(
1950  GetVertexMap(i), edgemodearray[j],
1951  -R(GetVertexMap(i), edgemodearray[j]));
1952  }
1953  for (j = 0; j < nfacemodestotal; ++j)
1954  {
1955  InvR.SetValue(
1956  GetVertexMap(i), facemodearray[j],
1957  -R(GetVertexMap(i), facemodearray[j]));
1958  for (n = 0; n < nedgemodestotal; ++n)
1959  {
1960  MatrixValue = InvR.GetValue(GetVertexMap(i),
1961  facemodearray[j])
1962  + R(GetVertexMap(i), edgemodearray[n])
1963  * R(edgemodearray[n], facemodearray[j]);
1964  InvR.SetValue(GetVertexMap(i),
1965  facemodearray[j],
1966  MatrixValue);
1967  }
1968  }
1969  }
1970 
1971  // Edge-face contributions
1972  for (i = 0; i < nedgemodestotal; ++i)
1973  {
1974  for (j = 0; j < nfacemodestotal; ++j)
1975  {
1976  InvR.SetValue(
1977  edgemodearray[i], facemodearray[j],
1978  -R(edgemodearray[i], facemodearray[j]));
1979  }
1980  }
1981 
1982  for (i = 0; i < nCoeffs; ++i)
1983  {
1984  InvR.SetValue(i, i, 1.0);
1985  }
1986 
1987  return inversetransformationmatrix;
1988  }
Array< OneD, unsigned int > GetTraceInverseBoundaryMap(int fid, StdRegions::Orientation faceOrient=StdRegions::eNoOrientation, int P1=-1, int P2=-1)
Array< OneD, unsigned int > GetEdgeInverseBoundaryMap(int eid)
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
NekMatrix< NekDouble, StandardMatrixTag > DNekMat
Definition: NekTypeDefs.hpp:51
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69

References Nektar::eFULL, and Vmath::Vcopy().

◆ v_BuildTransformationMatrix()

DNekMatSharedPtr Nektar::LocalRegions::Expansion3D::v_BuildTransformationMatrix ( const DNekScalMatSharedPtr r_bnd,
const StdRegions::MatrixType  matrixType 
)
protectedvirtual

The matrix component of \(\mathbf{R}\) is given by

\[ \mathbf{R^{T}_{v}}= -\mathbf{S}^{-1}_{ef,ef}\mathbf{S}^{T}_{v,ef}\]

For every vertex mode we extract the submatrices from statically condensed matrix \(\mathbf{S}\) corresponding to the coupling between the attached edges and faces of a vertex ( \(\mathbf{S_{ef,ef}}\)). This matrix is then inverted and multiplied by the submatrix representing the coupling between a vertex and the attached edges and faces ( \(\mathbf{S_{v,ef}}\)).

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1452 of file Expansion3D.cpp.

1455  {
1456  int nVerts, nEdges;
1457  int eid, fid, vid, n, i;
1458 
1459  int nBndCoeffs = NumBndryCoeffs();
1460 
1462 
1463  // Get geometric information about this element
1464  nVerts = GetNverts();
1465  nEdges = GetNedges();
1466 
1467  /*************************************/
1468  /* Vetex-edge & vertex-face matrices */
1469  /*************************************/
1470 
1471  /**
1472  * The matrix component of \f$\mathbf{R}\f$ is given by \f[
1473  * \mathbf{R^{T}_{v}}=
1474  * -\mathbf{S}^{-1}_{ef,ef}\mathbf{S}^{T}_{v,ef}\f]
1475  *
1476  * For every vertex mode we extract the submatrices from statically
1477  * condensed matrix \f$\mathbf{S}\f$ corresponding to the coupling
1478  * between the attached edges and faces of a vertex
1479  * (\f$\mathbf{S_{ef,ef}}\f$). This matrix is then inverted and
1480  * multiplied by the submatrix representing the coupling between a
1481  * vertex and the attached edges and faces
1482  * (\f$\mathbf{S_{v,ef}}\f$).
1483  */
1484 
1485  int nmodes;
1486  int m;
1487  NekDouble VertexEdgeFaceValue;
1488 
1489  // The number of connected edges/faces is 3 (for all elements)
1490  int nConnectedEdges = 3;
1491  int nConnectedFaces = 3;
1492 
1493  // Location in the matrix
1494  Array<OneD, Array<OneD, unsigned int> >
1495  MatEdgeLocation(nConnectedEdges);
1496  Array<OneD, Array<OneD, unsigned int> >
1497  MatFaceLocation(nConnectedFaces);
1498 
1499  // Define storage for vertex transpose matrix and zero all entries
1500  MatrixStorage storage = eFULL;
1501  DNekMatSharedPtr transformationmatrix;
1502 
1503  transformationmatrix =
1505  nBndCoeffs, nBndCoeffs, 0.0, storage);
1506 
1507  DNekMat &R = (*transformationmatrix);
1508 
1509 
1510  // Build the vertex-edge/face transform matrix: This matrix is
1511  // constructed from the submatrices corresponding to the couping
1512  // between each vertex and the attached edges/faces
1513  for (vid = 0; vid < nVerts; ++vid)
1514  {
1515  // Row and column size of the vertex-edge/face matrix
1516  int efRow =
1517  GetEdgeNcoeffs (geom->GetVertexEdgeMap(vid, 0)) +
1518  GetEdgeNcoeffs (geom->GetVertexEdgeMap(vid, 1)) +
1519  GetEdgeNcoeffs (geom->GetVertexEdgeMap(vid, 2)) +
1520  GetTraceIntNcoeffs(geom->GetVertexFaceMap(vid, 0)) +
1521  GetTraceIntNcoeffs(geom->GetVertexFaceMap(vid, 1)) +
1522  GetTraceIntNcoeffs(geom->GetVertexFaceMap(vid, 2)) - 6;
1523 
1524  int nedgemodesconnected =
1525  GetEdgeNcoeffs (geom->GetVertexEdgeMap(vid, 0)) +
1526  GetEdgeNcoeffs (geom->GetVertexEdgeMap(vid, 1)) +
1527  GetEdgeNcoeffs (geom->GetVertexEdgeMap(vid, 2)) - 6;
1528  Array<OneD, unsigned int> edgemodearray(nedgemodesconnected);
1529 
1530  int nfacemodesconnected =
1531  GetTraceIntNcoeffs(geom->GetVertexFaceMap(vid, 0)) +
1532  GetTraceIntNcoeffs(geom->GetVertexFaceMap(vid, 1)) +
1533  GetTraceIntNcoeffs(geom->GetVertexFaceMap(vid, 2));
1534  Array<OneD, unsigned int> facemodearray(nfacemodesconnected);
1535 
1536  int offset = 0;
1537 
1538  // Create array of edge modes
1539  for (eid = 0; eid < nConnectedEdges; ++eid)
1540  {
1541  MatEdgeLocation[eid] = GetEdgeInverseBoundaryMap(
1542  geom->GetVertexEdgeMap(vid, eid));
1543  nmodes = MatEdgeLocation[eid].size();
1544 
1545  if (nmodes)
1546  {
1547  Vmath::Vcopy(nmodes, &MatEdgeLocation[eid][0],
1548  1, &edgemodearray[offset], 1);
1549  }
1550 
1551  offset += nmodes;
1552  }
1553 
1554  offset = 0;
1555 
1556  // Create array of face modes
1557  for (fid = 0; fid < nConnectedFaces; ++fid)
1558  {
1559  MatFaceLocation[fid] = GetTraceInverseBoundaryMap(
1560  geom->GetVertexFaceMap(vid, fid));
1561  nmodes = MatFaceLocation[fid].size();
1562 
1563  if (nmodes)
1564  {
1565  Vmath::Vcopy(nmodes, &MatFaceLocation[fid][0],
1566  1, &facemodearray[offset], 1);
1567  }
1568  offset += nmodes;
1569  }
1570 
1571  DNekMatSharedPtr vertexedgefacetransformmatrix =
1573  1, efRow, 0.0, storage);
1574  DNekMat &Sveft = (*vertexedgefacetransformmatrix);
1575 
1576  DNekMatSharedPtr vertexedgefacecoupling =
1578  1, efRow, 0.0, storage);
1579  DNekMat &Svef = (*vertexedgefacecoupling);
1580 
1581  // Vertex-edge coupling
1582  for (n = 0; n < nedgemodesconnected; ++n)
1583  {
1584  // Matrix value for each coefficient location
1585  VertexEdgeFaceValue = (*r_bnd)(GetVertexMap(vid),
1586  edgemodearray[n]);
1587 
1588  // Set the value in the vertex edge/face matrix
1589  Svef.SetValue(0, n, VertexEdgeFaceValue);
1590  }
1591 
1592  // Vertex-face coupling
1593  for (n = 0; n < nfacemodesconnected; ++n)
1594  {
1595  // Matrix value for each coefficient location
1596  VertexEdgeFaceValue = (*r_bnd)(GetVertexMap(vid),
1597  facemodearray[n]);
1598 
1599  // Set the value in the vertex edge/face matrix
1600  Svef.SetValue(0, n + nedgemodesconnected,
1601  VertexEdgeFaceValue);
1602  }
1603 
1604  /*
1605  * Build the edge-face transform matrix: This matrix is
1606  * constructed from the submatrices corresponding to the couping
1607  * between the edges and faces on the attached faces/edges of a
1608  * vertex.
1609  */
1610 
1611  // Allocation of matrix to store edge/face-edge/face coupling
1612  DNekMatSharedPtr edgefacecoupling =
1614  efRow, efRow, 0.0, storage);
1615  DNekMat &Sefef = (*edgefacecoupling);
1616 
1617  NekDouble EdgeEdgeValue, FaceFaceValue;
1618 
1619  // Edge-edge coupling (S_{ee})
1620  for (m = 0; m < nedgemodesconnected; ++m)
1621  {
1622  for (n = 0; n < nedgemodesconnected; ++n)
1623  {
1624  // Matrix value for each coefficient location
1625  EdgeEdgeValue = (*r_bnd)(edgemodearray[n],
1626  edgemodearray[m]);
1627 
1628  // Set the value in the vertex edge/face matrix
1629  Sefef.SetValue(n, m, EdgeEdgeValue);
1630  }
1631  }
1632 
1633  // Face-face coupling (S_{ff})
1634  for (n = 0; n < nfacemodesconnected; ++n)
1635  {
1636  for (m = 0; m < nfacemodesconnected; ++m)
1637  {
1638  // Matrix value for each coefficient location
1639  FaceFaceValue = (*r_bnd)(facemodearray[n],
1640  facemodearray[m]);
1641  // Set the value in the vertex edge/face matrix
1642  Sefef.SetValue(nedgemodesconnected + n,
1643  nedgemodesconnected + m, FaceFaceValue);
1644  }
1645  }
1646 
1647  // Edge-face coupling (S_{ef} and trans(S_{ef}))
1648  for (n = 0; n < nedgemodesconnected; ++n)
1649  {
1650  for (m = 0; m < nfacemodesconnected; ++m)
1651  {
1652  // Matrix value for each coefficient location
1653  FaceFaceValue = (*r_bnd)(edgemodearray[n],
1654  facemodearray[m]);
1655 
1656  // Set the value in the vertex edge/face matrix
1657  Sefef.SetValue(n,
1658  nedgemodesconnected + m,
1659  FaceFaceValue);
1660 
1661  FaceFaceValue = (*r_bnd)(facemodearray[m],
1662  edgemodearray[n]);
1663 
1664  // and transpose
1665  Sefef.SetValue(nedgemodesconnected + m,
1666  n,
1667  FaceFaceValue);
1668  }
1669  }
1670 
1671  // Invert edge-face coupling matrix
1672  if (efRow)
1673  {
1674  Sefef.Invert();
1675 
1676 
1677  //R_{v}=-S_{v,ef}inv(S_{ef,ef})
1678  Sveft = -Svef * Sefef;
1679  }
1680 
1681  // Populate R with R_{ve} components
1682  for (n = 0; n < edgemodearray.size(); ++n)
1683  {
1684  R.SetValue(GetVertexMap(vid), edgemodearray[n],
1685  Sveft(0, n));
1686  }
1687 
1688  // Populate R with R_{vf} components
1689  for (n = 0; n < facemodearray.size(); ++n)
1690  {
1691  R.SetValue(GetVertexMap(vid), facemodearray[n],
1692  Sveft(0, n + nedgemodesconnected));
1693  }
1694  }
1695 
1696  /********************/
1697  /* edge-face matrix */
1698  /********************/
1699 
1700  /*
1701  * The matrix component of \f$\mathbf{R}\f$ is given by \f[
1702  * \mathbf{R^{T}_{ef}}=-\mathbf{S}^{-1}_{ff}\mathbf{S}^{T}_{ef}\f]
1703  *
1704  * For each edge extract the submatrices from statically condensed
1705  * matrix \f$\mathbf{S}\f$ corresponding to inner products of modes
1706  * on the two attached faces within themselves as well as the
1707  * coupling matrix between the two faces
1708  * (\f$\mathbf{S}_{ff}\f$). This matrix of face coupling is then
1709  * inverted and multiplied by the submatrices of corresponding to
1710  * the coupling between the edge and attached faces
1711  * (\f$\mathbf{S}_{ef}\f$).
1712  */
1713 
1714  NekDouble EdgeFaceValue, FaceFaceValue;
1715  int efCol, efRow, nedgemodes;
1716 
1717  // Number of attached faces is always 2
1718  nConnectedFaces = 2;
1719 
1720  // Location in the matrix
1721  MatEdgeLocation = Array<OneD, Array<OneD, unsigned int> >
1722  (nEdges);
1723  MatFaceLocation = Array<OneD, Array<OneD, unsigned int> >
1724  (nConnectedFaces);
1725 
1726  // Build the edge/face transform matrix: This matrix is constructed
1727  // from the submatrices corresponding to the couping between a
1728  // specific edge and the two attached faces.
1729  for (eid = 0; eid < nEdges; ++eid)
1730  {
1731  // Row and column size of the vertex-edge/face matrix
1732  efCol = GetTraceIntNcoeffs(geom->GetEdgeFaceMap(eid, 0)) +
1733  GetTraceIntNcoeffs(geom->GetEdgeFaceMap(eid, 1));
1734  efRow = GetEdgeNcoeffs(eid) - 2;
1735 
1736  // Edge-face coupling matrix
1737  DNekMatSharedPtr efedgefacecoupling =
1739  efRow, efCol, 0.0, storage);
1740  DNekMat &Mef = (*efedgefacecoupling);
1741 
1742  // Face-face coupling matrix
1743  DNekMatSharedPtr effacefacecoupling =
1745  efCol, efCol, 0.0, storage);
1746  DNekMat &Meff = (*effacefacecoupling);
1747 
1748  // Edge-face transformation matrix
1749  DNekMatSharedPtr edgefacetransformmatrix =
1751  efRow, efCol, 0.0, storage);
1752  DNekMat &Meft = (*edgefacetransformmatrix);
1753 
1754  int nfacemodesconnected =
1755  GetTraceIntNcoeffs(geom->GetEdgeFaceMap(eid, 0)) +
1756  GetTraceIntNcoeffs(geom->GetEdgeFaceMap(eid, 1));
1757  Array<OneD, unsigned int>
1758  facemodearray(nfacemodesconnected);
1759 
1760  // Create array of edge modes
1761  Array<OneD, unsigned int> inedgearray
1763  nedgemodes = GetEdgeNcoeffs(eid) - 2;
1764  Array<OneD, unsigned int> edgemodearray(nedgemodes);
1765 
1766  if (nedgemodes)
1767  {
1768  Vmath::Vcopy(nedgemodes, &inedgearray[0],
1769  1, &edgemodearray[0], 1);
1770  }
1771 
1772  int offset = 0;
1773 
1774  // Create array of face modes
1775  for (fid = 0; fid < nConnectedFaces; ++fid)
1776  {
1777  MatFaceLocation[fid] = GetTraceInverseBoundaryMap(
1778  geom->GetEdgeFaceMap(eid, fid));
1779  nmodes = MatFaceLocation[fid].size();
1780 
1781  if (nmodes)
1782  {
1783  Vmath::Vcopy(nmodes, &MatFaceLocation[fid][0],
1784  1, &facemodearray[offset], 1);
1785  }
1786  offset += nmodes;
1787  }
1788 
1789  // Edge-face coupling
1790  for (n = 0; n < nedgemodes; ++n)
1791  {
1792  for (m = 0; m < nfacemodesconnected; ++m)
1793  {
1794  // Matrix value for each coefficient location
1795  EdgeFaceValue = (*r_bnd)(edgemodearray[n],
1796  facemodearray[m]);
1797 
1798  // Set the value in the edge/face matrix
1799  Mef.SetValue(n, m, EdgeFaceValue);
1800  }
1801  }
1802 
1803  // Face-face coupling
1804  for (n = 0; n < nfacemodesconnected; ++n)
1805  {
1806  for (m = 0; m < nfacemodesconnected; ++m)
1807  {
1808  // Matrix value for each coefficient location
1809  FaceFaceValue = (*r_bnd)(facemodearray[n],
1810  facemodearray[m]);
1811 
1812  // Set the value in the vertex edge/face matrix
1813  Meff.SetValue(n, m, FaceFaceValue);
1814  }
1815  }
1816 
1817  if (efCol)
1818  {
1819  // Invert edge-face coupling matrix
1820  Meff.Invert();
1821 
1822  // trans(R_{ef})=-S_{ef}*(inv(S_{ff})
1823  Meft = -Mef * Meff;
1824  }
1825 
1826  //Populate transformation matrix with Meft
1827  for (n = 0; n < Meft.GetRows(); ++n)
1828  {
1829  for (m = 0; m < Meft.GetColumns(); ++m)
1830  {
1831  R.SetValue(edgemodearray[n], facemodearray[m],
1832  Meft(n, m));
1833  }
1834  }
1835  }
1836 
1837  for (i = 0; i < R.GetRows(); ++i)
1838  {
1839  R.SetValue(i, i, 1.0);
1840  }
1841 
1842  if ((matrixType == StdRegions::ePreconR)||
1843  (matrixType == StdRegions::ePreconRMass))
1844  {
1845  return transformationmatrix;
1846  }
1847  else
1848  {
1849  NEKERROR(ErrorUtil::efatal, "unkown matrix type" );
1850  return NullDNekMatSharedPtr;
1851  }
1852  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
static DNekMatSharedPtr NullDNekMatSharedPtr
Definition: NekTypeDefs.hpp:78

References Nektar::eFULL, Nektar::StdRegions::ePreconR, Nektar::StdRegions::ePreconRMass, NEKERROR, Nektar::NullDNekMatSharedPtr, and Vmath::Vcopy().

◆ v_BuildVertexMatrix()

DNekMatSharedPtr Nektar::LocalRegions::Expansion3D::v_BuildVertexMatrix ( const DNekScalMatSharedPtr r_bnd)
protectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1421 of file Expansion3D.cpp.

1423  {
1424  MatrixStorage storage = eFULL;
1425  DNekMatSharedPtr vertexmatrix;
1426 
1427  int nVerts, vid1, vid2, vMap1, vMap2;
1428  NekDouble VertexValue;
1429 
1430  nVerts = GetNverts();
1431 
1432  vertexmatrix =
1434  nVerts, nVerts, 0.0, storage);
1435  DNekMat &VertexMat = (*vertexmatrix);
1436 
1437  for (vid1 = 0; vid1 < nVerts; ++vid1)
1438  {
1439  vMap1 = GetVertexMap(vid1,true);
1440 
1441  for (vid2 = 0; vid2 < nVerts; ++vid2)
1442  {
1443  vMap2 = GetVertexMap(vid2,true);
1444  VertexValue = (*r_bnd)(vMap1, vMap2);
1445  VertexMat.SetValue(vid1, vid2, VertexValue);
1446  }
1447  }
1448 
1449  return vertexmatrix;
1450  }

References Nektar::eFULL.

◆ v_DGDeriv()

void Nektar::LocalRegions::Expansion3D::v_DGDeriv ( const int  dir,
const Array< OneD, const NekDouble > &  incoeffs,
Array< OneD, ExpansionSharedPtr > &  FaceExp,
Array< OneD, Array< OneD, NekDouble > > &  faceCoeffs,
Array< OneD, NekDouble > &  out_d 
)
protectedvirtual

Evaluate coefficients of weak deriviative in the direction dir given the input coefficicents incoeffs and the imposed boundary values in EdgeExp (which will have its phys space updated).

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1231 of file Expansion3D.cpp.

1237  {
1238  int ncoeffs = GetNcoeffs();
1242 
1244  DNekScalMat &Dmat = *GetLocMatrix(DerivType[dir]);
1245 
1246  Array<OneD, NekDouble> coeffs = incoeffs;
1247  DNekVec Coeffs (ncoeffs,coeffs, eWrapper);
1248 
1249  Coeffs = Transpose(Dmat)*Coeffs;
1250  Vmath::Neg(ncoeffs, coeffs,1);
1251 
1252  // Add the boundary integral including the relevant part of
1253  // the normal
1254  AddNormTraceInt(dir, FaceExp, faceCoeffs, coeffs);
1255 
1256  DNekVec Out_d (ncoeffs,out_d,eWrapper);
1257 
1258  Out_d = InvMass*Coeffs;
1259  }
void AddNormTraceInt(const int dir, Array< OneD, ExpansionSharedPtr > &FaceExp, Array< OneD, Array< OneD, NekDouble > > &faceCoeffs, Array< OneD, NekDouble > &outarray)
NekMatrix< InnerMatrixType, BlockMatrixTag > Transpose(NekMatrix< InnerMatrixType, BlockMatrixTag > &rhs)
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:461

References Nektar::StdRegions::eInvMass, Nektar::StdRegions::eWeakDeriv0, Nektar::StdRegions::eWeakDeriv1, Nektar::StdRegions::eWeakDeriv2, Nektar::eWrapper, Vmath::Neg(), and Nektar::Transpose().

◆ v_GenMatrix()

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

Computes matrices needed for the HDG formulation. References to equations relate to the following paper (with a suitable changes in formulation to adapt to 3D): 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 NOTE: VARIABLE COEFFICIENTS CASE IS NOT IMPLEMENTED

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TetExp, Nektar::LocalRegions::PyrExp, Nektar::LocalRegions::PrismExp, and Nektar::LocalRegions::HexExp.

Definition at line 451 of file Expansion3D.cpp.

452  {
453  //Variable coefficients are not implemented/////////
454  ASSERTL1(!mkey.HasVarCoeff(StdRegions::eVarCoeffD00),
455  "Matrix construction is not implemented for variable "
456  "coefficients at the moment");
457  ////////////////////////////////////////////////////
458 
459  DNekMatSharedPtr returnval;
460 
461  switch(mkey.GetMatrixType())
462  {
463  // (Z^e)^{-1} (Eqn. 33, P22)
465  {
467  "HybridDGHelmholtz matrix not set up "
468  "for non boundary-interior expansions");
469 
470  int i,j,k;
471  NekDouble lambdaval = mkey.GetConstFactor(StdRegions::eFactorLambda);
472  NekDouble tau = mkey.GetConstFactor(StdRegions::eFactorTau);
473  int ncoeffs = GetNcoeffs();
474  int nfaces = GetNtraces();
475 
476  Array<OneD,unsigned int> fmap;
477  Array<OneD,int> sign;
478  ExpansionSharedPtr FaceExp;
479  ExpansionSharedPtr FaceExp2;
480 
481  int order_f, coordim = GetCoordim();
486 
487  returnval = MemoryManager<DNekMat>::AllocateSharedPtr(ncoeffs,ncoeffs);
488  DNekMat &Mat = *returnval;
489  Vmath::Zero(ncoeffs*ncoeffs,Mat.GetPtr(),1);
490 
491  // StdRegions::VarCoeffType Coeffs[3] = {StdRegions::eVarCoeffD00,
492  // StdRegions::eVarCoeffD11,
493  // StdRegions::eVarCoeffD22};
494  StdRegions::VarCoeffMap::const_iterator x;
495  const StdRegions::VarCoeffMap &varcoeffs =
496  mkey.GetVarCoeffs();
497 
498  for(i = 0; i < coordim; ++i)
499  {
500  if ((x = varcoeffs.find(StdRegions::eVarCoeffMF1x))
501  != varcoeffs.end())
502  {
503  StdRegions::VarCoeffMap VarCoeffDirDeriv;
504  VarCoeffDirDeriv[StdRegions::eVarCoeffMF] =
505  v_GetMF(i,coordim,varcoeffs);
506  VarCoeffDirDeriv[StdRegions::eVarCoeffMFDiv] =
507  v_GetMFDiv(i,varcoeffs);
508 
509  MatrixKey Dmatkey(StdRegions::eWeakDirectionalDeriv,
510  DetShapeType(), *this,
512  VarCoeffDirDeriv);
513 
514  DNekScalMat &Dmat = *GetLocMatrix(Dmatkey);
515 
517  Weight[StdRegions::eVarCoeffMass] =
518  v_GetMFMag(i,mkey.GetVarCoeffs());
519 
520  MatrixKey invMasskey(StdRegions::eInvMass,
521  DetShapeType(), *this,
523  Weight);
524 
525  DNekScalMat &invMass = *GetLocMatrix(invMasskey);
526 
527  Mat = Mat + Dmat*invMass*Transpose(Dmat);
528  }
529  else
530  {
531  DNekScalMat &Dmat = *GetLocMatrix(DerivType[i]);
532  Mat = Mat + Dmat*invMass*Transpose(Dmat);
533  }
534 
535  /*
536  if(mkey.HasVarCoeff(Coeffs[i]))
537  {
538  MatrixKey DmatkeyL(DerivType[i], DetExpansionType(), *this,
539  StdRegions::NullConstFactorMap,
540  mkey.GetVarCoeffAsMap(Coeffs[i]));
541  MatrixKey DmatkeyR(DerivType[i], DetExpansionType(), *this);
542 
543  DNekScalMat &DmatL = *GetLocMatrix(DmatkeyL);
544  DNekScalMat &DmatR = *GetLocMatrix(DmatkeyR);
545  Mat = Mat + DmatL*invMass*Transpose(DmatR);
546  }
547  else
548  {
549  DNekScalMat &Dmat = *GetLocMatrix(DerivType[i]);
550  Mat = Mat + Dmat*invMass*Transpose(Dmat);
551  }
552  */
553  }
554 
555  // Add Mass Matrix Contribution for Helmholtz problem
557  Mat = Mat + lambdaval*Mass;
558 
559  // Add tau*E_l using elemental mass matrices on each edge
560  for(i = 0; i < nfaces; ++i)
561  {
562  FaceExp = GetTraceExp(i);
563  order_f = FaceExp->GetNcoeffs();
564 
565  IndexMapKey ikey(eFaceToElement, DetShapeType(),
567  GetBasisNumModes(2), i, GetTraceOrient(i));
569 
570  // @TODO: Document
571  /*
572  StdRegions::VarCoeffMap edgeVarCoeffs;
573  if (mkey.HasVarCoeff(StdRegions::eVarCoeffD00))
574  {
575  Array<OneD, NekDouble> mu(nq);
576  GetPhysEdgeVarCoeffsFromElement(
577  i, EdgeExp2,
578  mkey.GetVarCoeff(StdRegions::eVarCoeffD00), mu);
579  edgeVarCoeffs[StdRegions::eVarCoeffMass] = mu;
580  }
581  DNekScalMat &eMass = *EdgeExp->GetLocMatrix(
582  StdRegions::eMass,
583  StdRegions::NullConstFactorMap, edgeVarCoeffs);
584  */
585 
586  DNekScalMat &eMass = *FaceExp->GetLocMatrix(StdRegions::eMass);
587 
588  for(j = 0; j < order_f; ++j)
589  {
590  for(k = 0; k < order_f; ++k)
591  {
592  Mat((*map)[j].index,(*map)[k].index) +=
593  tau*(*map)[j].sign*(*map)[k].sign*eMass(j,k);
594  }
595  }
596  }
597  break;
598  }
599 
600  // U^e (P22)
602  {
603  int i,j,k;
604  int nbndry = NumDGBndryCoeffs();
605  int ncoeffs = GetNcoeffs();
606  int nfaces = GetNtraces();
607  NekDouble tau = mkey.GetConstFactor(StdRegions::eFactorTau);
608 
609  Array<OneD,NekDouble> lambda(nbndry);
610  DNekVec Lambda(nbndry,lambda,eWrapper);
611  Array<OneD,NekDouble> ulam(ncoeffs);
612  DNekVec Ulam(ncoeffs,ulam,eWrapper);
613  Array<OneD,NekDouble> f(ncoeffs);
614  DNekVec F(ncoeffs,f,eWrapper);
615 
616  ExpansionSharedPtr FaceExp;
617  // declare matrix space
618  returnval = MemoryManager<DNekMat>::AllocateSharedPtr(ncoeffs,nbndry);
619  DNekMat &Umat = *returnval;
620 
621  // Z^e matrix
622  MatrixKey newkey(StdRegions::eInvHybridDGHelmholtz, DetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
623  DNekScalMat &invHmat = *GetLocMatrix(newkey);
624 
625  Array<OneD,unsigned int> fmap;
626  Array<OneD,int> sign;
627 
628  //alternative way to add boundary terms contribution
629  int bndry_cnt = 0;
630  for(i = 0; i < nfaces; ++i)
631  {
632  FaceExp = GetTraceExp(i);//temporary, need to rewrite AddHDGHelmholtzFaceTerms
633  int nface = GetTraceNcoeffs(i);
634  Array<OneD, NekDouble> face_lambda(nface);
635 
636  const Array<OneD, const Array<OneD, NekDouble> > normals
637  = GetTraceNormal(i);
638 
639  for(j = 0; j < nface; ++j)
640  {
641  Vmath::Zero(nface,&face_lambda[0],1);
642  Vmath::Zero(ncoeffs,&f[0],1);
643  face_lambda[j] = 1.0;
644 
645  SetFaceToGeomOrientation(i, face_lambda);
646 
647  Array<OneD, NekDouble> tmp(FaceExp->GetTotPoints());
648  FaceExp->BwdTrans(face_lambda, tmp);
649  AddHDGHelmholtzFaceTerms(tau, i, tmp, mkey.GetVarCoeffs(), f);
650 
651  Ulam = invHmat*F; // generate Ulam from lambda
652 
653  // fill column of matrix
654  for(k = 0; k < ncoeffs; ++k)
655  {
656  Umat(k,bndry_cnt) = Ulam[k];
657  }
658 
659  ++bndry_cnt;
660  }
661  }
662 
663  //// Set up face expansions from local geom info
664  //for(i = 0; i < nfaces; ++i)
665  //{
666  // FaceExp[i] = GetTraceExp(i);
667  //}
668  //
669  //// for each degree of freedom of the lambda space
670  //// calculate Umat entry
671  //// Generate Lambda to U_lambda matrix
672  //for(j = 0; j < nbndry; ++j)
673  //{
674  // // standard basis vectors e_j
675  // Vmath::Zero(nbndry,&lambda[0],1);
676  // Vmath::Zero(ncoeffs,&f[0],1);
677  // lambda[j] = 1.0;
678  //
679  // //cout << Lambda;
680  // SetTraceToGeomOrientation(lambda);
681  // //cout << Lambda << endl;
682  //
683  // // Compute F = [I D_1 M^{-1} D_2 M^{-1}] C e_j
684  // AddHDGHelmholtzTraceTerms(tau, lambda, FaceExp, mkey.GetVarCoeffs(), f);
685  //
686  // // Compute U^e_j
687  // Ulam = invHmat*F; // generate Ulam from lambda
688  //
689  // // fill column of matrix
690  // for(k = 0; k < ncoeffs; ++k)
691  // {
692  // Umat(k,j) = Ulam[k];
693  // }
694  //}
695  }
696  break;
697  // Q_0, Q_1, Q_2 matrices (P23)
698  // Each are a product of a row of Eqn 32 with the C matrix.
699  // Rather than explicitly computing all of Eqn 32, we note each
700  // row is almost a multiple of U^e, so use that as our starting
701  // point.
705  {
706  int i = 0;
707  int j = 0;
708  int k = 0;
709  int dir = 0;
710  int nbndry = NumDGBndryCoeffs();
711  int coordim = GetCoordim();
712  int ncoeffs = GetNcoeffs();
713  int nfaces = GetNtraces();
714 
715  Array<OneD,NekDouble> lambda(nbndry);
716  DNekVec Lambda(nbndry,lambda,eWrapper);
717  Array<OneD, ExpansionSharedPtr> FaceExp(nfaces);
718 
719  Array<OneD,NekDouble> ulam(ncoeffs);
720  DNekVec Ulam(ncoeffs,ulam,eWrapper);
721  Array<OneD,NekDouble> f(ncoeffs);
722  DNekVec F(ncoeffs,f,eWrapper);
723 
724  // declare matrix space
725  returnval = MemoryManager<DNekMat>::AllocateSharedPtr(ncoeffs,nbndry);
726  DNekMat &Qmat = *returnval;
727 
728  // Lambda to U matrix
729  MatrixKey lamToUkey(StdRegions::eHybridDGLamToU, DetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
730  DNekScalMat &lamToU = *GetLocMatrix(lamToUkey);
731 
732  // Inverse mass matrix
734 
735  for(i = 0; i < nfaces; ++i)
736  {
737  FaceExp[i] = GetTraceExp(i);
738  }
739 
740  //Weak Derivative matrix
742  switch(mkey.GetMatrixType())
743  {
745  dir = 0;
747  break;
749  dir = 1;
751  break;
753  dir = 2;
755  break;
756  default:
757  ASSERTL0(false,"Direction not known");
758  break;
759  }
760 
761  //DNekScalMatSharedPtr Dmat;
762  //DNekScalMatSharedPtr &invMass;
763  StdRegions::VarCoeffMap::const_iterator x;
764  const StdRegions::VarCoeffMap &varcoeffs =
765  mkey.GetVarCoeffs();
766  if ((x = varcoeffs.find(StdRegions::eVarCoeffMF1x)) !=
767  varcoeffs.end())
768  {
769  StdRegions::VarCoeffMap VarCoeffDirDeriv;
770  VarCoeffDirDeriv[StdRegions::eVarCoeffMF] =
771  v_GetMF(dir,coordim,varcoeffs);
772  VarCoeffDirDeriv[StdRegions::eVarCoeffMFDiv] =
773  v_GetMFDiv(dir,varcoeffs);
774 
775  MatrixKey Dmatkey(StdRegions::eWeakDirectionalDeriv,
776  DetShapeType(), *this,
778  VarCoeffDirDeriv);
779 
780  Dmat = GetLocMatrix(Dmatkey);
781 
783  Weight[StdRegions::eVarCoeffMass] =
784  v_GetMFMag(dir,mkey.GetVarCoeffs());
785 
786  MatrixKey invMasskey(StdRegions::eInvMass,
787  DetShapeType(), *this,
789  Weight);
790 
791  invMass = *GetLocMatrix(invMasskey);
792  }
793  else
794  {
795  // Inverse mass matrix
797  }
798 
799  // for each degree of freedom of the lambda space
800  // calculate Qmat entry
801  // Generate Lambda to Q_lambda matrix
802  for(j = 0; j < nbndry; ++j)
803  {
804  Vmath::Zero(nbndry,&lambda[0],1);
805  lambda[j] = 1.0;
806 
807  // for lambda[j] = 1 this is the solution to ulam
808  for(k = 0; k < ncoeffs; ++k)
809  {
810  Ulam[k] = lamToU(k,j);
811  }
812 
813  // -D^T ulam
814  Vmath::Neg(ncoeffs,&ulam[0],1);
815  F = Transpose(*Dmat)*Ulam;
816 
818 
819  // Add the C terms resulting from the I's on the
820  // diagonals of Eqn 32
821  AddNormTraceInt(dir,lambda,FaceExp,f,mkey.GetVarCoeffs());
822 
823  // finally multiply by inverse mass matrix
824  Ulam = invMass*F;
825 
826  // fill column of matrix (Qmat is in column major format)
827  Vmath::Vcopy(ncoeffs,&ulam[0],1,&(Qmat.GetPtr())[0]+j*ncoeffs,1);
828  }
829  }
830  break;
831  // Matrix K (P23)
833  {
834  int i,j,f,cnt;
835  int order_f, nquad_f;
836  int nbndry = NumDGBndryCoeffs();
837  int nfaces = GetNtraces();
838  NekDouble tau = mkey.GetConstFactor(StdRegions::eFactorTau);
839 
840  Array<OneD,NekDouble> work, varcoeff_work;
841  Array<OneD,const Array<OneD, NekDouble> > normals;
842  Array<OneD, ExpansionSharedPtr> FaceExp(nfaces);
843  Array<OneD, NekDouble> lam(nbndry);
844 
845  Array<OneD,unsigned int> fmap;
846  Array<OneD, int> sign;
847 
848  // declare matrix space
849  returnval = MemoryManager<DNekMat>::AllocateSharedPtr(nbndry, nbndry);
850  DNekMat &BndMat = *returnval;
851 
852  DNekScalMatSharedPtr LamToQ[3];
853 
854  // Matrix to map Lambda to U
855  MatrixKey LamToUkey(StdRegions::eHybridDGLamToU, DetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
856  DNekScalMat &LamToU = *GetLocMatrix(LamToUkey);
857 
858  // Matrix to map Lambda to Q0
859  MatrixKey LamToQ0key(StdRegions::eHybridDGLamToQ0, DetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
860  LamToQ[0] = GetLocMatrix(LamToQ0key);
861 
862  // Matrix to map Lambda to Q1
863  MatrixKey LamToQ1key(StdRegions::eHybridDGLamToQ1, DetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
864  LamToQ[1] = GetLocMatrix(LamToQ1key);
865 
866  // Matrix to map Lambda to Q2
867  MatrixKey LamToQ2key(StdRegions::eHybridDGLamToQ2, DetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
868  LamToQ[2] = GetLocMatrix(LamToQ2key);
869 
870  // Set up edge segment expansions from local geom info
871  const StdRegions::VarCoeffMap &varcoeffs = mkey.GetVarCoeffs();
872  for(i = 0; i < nfaces; ++i)
873  {
874  FaceExp[i] = GetTraceExp(i);
875  }
876 
877  // Set up matrix derived from <mu, Q_lam.n - \tau (U_lam - Lam) >
878  for(i = 0; i < nbndry; ++i)
879  {
880  cnt = 0;
881 
882  Vmath::Zero(nbndry,lam,1);
883  lam[i] = 1.0;
885 
886  for(f = 0; f < nfaces; ++f)
887  {
888  order_f = FaceExp[f]->GetNcoeffs();
889  nquad_f = FaceExp[f]->GetNumPoints(0)*FaceExp[f]->GetNumPoints(1);
890  normals = GetTraceNormal(f);
891 
892  work = Array<OneD,NekDouble>(nquad_f);
893  varcoeff_work = Array<OneD, NekDouble>(nquad_f);
894 
895  IndexMapKey ikey(eFaceToElement, DetShapeType(),
897  GetBasisNumModes(2), f, GetTraceOrient(f));
899 
900  // @TODO Variable coefficients
901  /*
902  StdRegions::VarCoeffType VarCoeff[3] = {StdRegions::eVarCoeffD00,
903  StdRegions::eVarCoeffD11,
904  StdRegions::eVarCoeffD22};
905  const StdRegions::VarCoeffMap &varcoeffs = mkey.GetVarCoeffs();
906  StdRegions::VarCoeffMap::const_iterator x;
907  */
908 
909  // Q0 * n0 (BQ_0 terms)
910  Array<OneD, NekDouble> faceCoeffs(order_f);
911  Array<OneD, NekDouble> facePhys (nquad_f);
912  for(j = 0; j < order_f; ++j)
913  {
914  faceCoeffs[j] = (*map)[j].sign*(*LamToQ[0])((*map)[j].index,i);
915  }
916 
917  FaceExp[f]->BwdTrans(faceCoeffs, facePhys);
918 
919  // @TODO Variable coefficients
920  // Multiply by variable coefficient
921  /*
922  if ((x = varcoeffs.find(VarCoeff[0])) != varcoeffs.end())
923  {
924  GetPhysEdgeVarCoeffsFromElement(e,EdgeExp[e],x->second,varcoeff_work);
925  Vmath::Vmul(nquad_e,varcoeff_work,1,EdgeExp[e]->GetPhys(),1,EdgeExp[e]->UpdatePhys(),1);
926  }
927  */
928 
929  if (varcoeffs.find(StdRegions::eVarCoeffMF1x)
930  != varcoeffs.end())
931  {
932  Array<OneD, NekDouble> ncdotMF =
933  v_GetnFacecdotMF(0, f, FaceExp[f],
934  normals, varcoeffs);
935 
936  Vmath::Vmul(nquad_f, ncdotMF, 1,
937  facePhys, 1,
938  work, 1);
939  }
940  else
941  {
942  Vmath::Vmul(nquad_f, normals[0], 1,
943  facePhys, 1,
944  work, 1);
945  }
946 
947  // Q1 * n1 (BQ_1 terms)
948  for(j = 0; j < order_f; ++j)
949  {
950  faceCoeffs[j] = (*map)[j].sign*(*LamToQ[1])((*map)[j].index,i);
951  }
952 
953  FaceExp[f]->BwdTrans(faceCoeffs, facePhys);
954 
955  // @TODO Variable coefficients
956  // Multiply by variable coefficients
957  /*
958  if ((x = varcoeffs.find(VarCoeff[1])) != varcoeffs.end())
959  {
960  GetPhysEdgeVarCoeffsFromElement(e,EdgeExp[e],x->second,varcoeff_work);
961  Vmath::Vmul(nquad_e,varcoeff_work,1,EdgeExp[e]->GetPhys(),1,EdgeExp[e]->UpdatePhys(),1);
962  }
963  */
964 
965  if ((varcoeffs.find(StdRegions::eVarCoeffMF1x)) !=
966  varcoeffs.end())
967  {
968  Array<OneD, NekDouble> ncdotMF =
969  v_GetnFacecdotMF(1, f, FaceExp[f],
970  normals, varcoeffs);
971 
972  Vmath::Vvtvp(nquad_f, ncdotMF, 1,
973  facePhys, 1,
974  work, 1,
975  work, 1);
976  }
977  else
978  {
979  Vmath::Vvtvp(nquad_f, normals[1], 1,
980  facePhys, 1,
981  work, 1,
982  work, 1);
983  }
984 
985  // Q2 * n2 (BQ_2 terms)
986  for(j = 0; j < order_f; ++j)
987  {
988  faceCoeffs[j] = (*map)[j].sign*(*LamToQ[2])((*map)[j].index,i);
989  }
990 
991  FaceExp[f]->BwdTrans(faceCoeffs, facePhys);
992 
993  // @TODO Variable coefficients
994  // Multiply by variable coefficients
995  /*
996  if ((x = varcoeffs.find(VarCoeff[2])) != varcoeffs.end())
997  {
998  GetPhysEdgeVarCoeffsFromElement(e,EdgeExp[e],x->second,varcoeff_work);
999  Vmath::Vmul(nquad_e,varcoeff_work,1,EdgeExp[e]->GetPhys(),1,EdgeExp[e]->UpdatePhys(),1);
1000  }
1001  */
1002 
1003  if (varcoeffs.find(StdRegions::eVarCoeffMF1x) !=
1004  varcoeffs.end())
1005  {
1006  Array<OneD, NekDouble> ncdotMF =
1007  v_GetnFacecdotMF(2, f, FaceExp[f],
1008  normals, varcoeffs);
1009 
1010  Vmath::Vvtvp(nquad_f, ncdotMF, 1,
1011  facePhys, 1,
1012  work, 1,
1013  work, 1);
1014  }
1015  else
1016  {
1017  Vmath::Vvtvp(nquad_f, normals[2], 1,
1018  facePhys, 1,
1019  work, 1,
1020  work, 1);
1021  }
1022 
1023  // - tau (ulam - lam)
1024  // Corresponds to the G and BU terms.
1025  for(j = 0; j < order_f; ++j)
1026  {
1027  faceCoeffs[j] = (*map)[j].sign*LamToU((*map)[j].index,i) - lam[cnt+j];
1028  }
1029 
1030  FaceExp[f]->BwdTrans(faceCoeffs, facePhys);
1031 
1032  // @TODO Variable coefficients
1033  // Multiply by variable coefficients
1034  /*
1035  if ((x = varcoeffs.find(VarCoeff[0])) != varcoeffs.end())
1036  {
1037  GetPhysEdgeVarCoeffsFromElement(e,FaceExp[f],x->second,varcoeff_work);
1038  Vmath::Vmul(nquad_f,varcoeff_work,1,FaceExp[f]->GetPhys(),1,FaceExp[f]->UpdatePhys(),1);
1039  }
1040  */
1041 
1042  Vmath::Svtvp(nquad_f, -tau, facePhys, 1,
1043  work, 1, work, 1);
1044 
1045  // @TODO Add variable coefficients
1046  FaceExp[f]->IProductWRTBase(work, faceCoeffs);
1047 
1048  SetFaceToGeomOrientation(f, faceCoeffs);
1049 
1050  for(j = 0; j < order_f; ++j)
1051  {
1052  BndMat(cnt+j,i) = faceCoeffs[j];
1053  }
1054 
1055  cnt += order_f;
1056  }
1057  }
1058  }
1059  break;
1060  //HDG postprocessing
1062  {
1063  MatrixKey lapkey(StdRegions::eLaplacian, DetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
1064  DNekScalMat &LapMat = *GetLocMatrix(lapkey);
1065 
1066  returnval = MemoryManager<DNekMat>::AllocateSharedPtr(LapMat.GetRows(),LapMat.GetColumns());
1067  DNekMatSharedPtr lmat = returnval;
1068 
1069  (*lmat) = LapMat;
1070 
1071  // replace first column with inner product wrt 1
1072  int nq = GetTotPoints();
1073  Array<OneD, NekDouble> tmp(nq);
1074  Array<OneD, NekDouble> outarray(m_ncoeffs);
1075  Vmath::Fill(nq,1.0,tmp,1);
1076  IProductWRTBase(tmp, outarray);
1077 
1078  Vmath::Vcopy(m_ncoeffs,&outarray[0],1,
1079  &(lmat->GetPtr())[0],1);
1080 
1081  lmat->Invert();
1082  }
1083  break;
1084  default:
1085  ASSERTL0(false,"This matrix type cannot be generated from this class");
1086  break;
1087  }
1088 
1089  return returnval;
1090  }
void SetTraceToGeomOrientation(Array< OneD, NekDouble > &inout)
Align trace orientation with the geometry orientation.
void AddHDGHelmholtzFaceTerms(const NekDouble tau, const int edge, Array< OneD, NekDouble > &facePhys, const StdRegions::VarCoeffMap &dirForcing, Array< OneD, NekDouble > &outarray)
Definition: Expansion3D.cpp:52
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:134
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:537
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
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:565
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:513
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:436
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45

References ASSERTL0, ASSERTL1, Nektar::LocalRegions::eFaceToElement, Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eFactorTau, Nektar::StdRegions::eHybridDGHelmBndLam, Nektar::StdRegions::eHybridDGHelmholtz, Nektar::StdRegions::eHybridDGLamToQ0, Nektar::StdRegions::eHybridDGLamToQ1, Nektar::StdRegions::eHybridDGLamToQ2, Nektar::StdRegions::eHybridDGLamToU, Nektar::StdRegions::eInvHybridDGHelmholtz, Nektar::StdRegions::eInvLaplacianWithUnityMean, Nektar::StdRegions::eInvMass, Nektar::StdRegions::eLaplacian, Nektar::StdRegions::eMass, Nektar::StdRegions::eVarCoeffD00, 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::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdMatrixKey::GetConstFactors(), Nektar::StdRegions::StdMatrixKey::GetMatrixType(), Nektar::StdRegions::StdMatrixKey::GetVarCoeffs(), Nektar::StdRegions::StdMatrixKey::HasVarCoeff(), Vmath::Neg(), Nektar::StdRegions::NullConstFactorMap, sign, Vmath::Svtvp(), Nektar::Transpose(), Vmath::Vcopy(), Vmath::Vmul(), Vmath::Vvtvp(), and Vmath::Zero().

Referenced by Nektar::LocalRegions::HexExp::v_GenMatrix(), Nektar::LocalRegions::PrismExp::v_GenMatrix(), Nektar::LocalRegions::PyrExp::v_GenMatrix(), and Nektar::LocalRegions::TetExp::v_GenMatrix().

◆ v_GetnFacecdotMF()

Array< OneD, NekDouble > Nektar::LocalRegions::Expansion3D::v_GetnFacecdotMF ( const int  dir,
const int  face,
ExpansionSharedPtr FaceExp_f,
const Array< OneD, const Array< OneD, NekDouble > > &  normals,
const StdRegions::VarCoeffMap varcoeffs 
)
protectedvirtual

Definition at line 2443 of file Expansion3D.cpp.

2449  {
2450  int nquad_f = FaceExp_f->GetNumPoints(0)*FaceExp_f->GetNumPoints(1);
2451  int coordim = GetCoordim();
2452 
2453  int nquad0 = m_base[0]->GetNumPoints();
2454  int nquad1 = m_base[1]->GetNumPoints();
2455  int nquad2 = m_base[2]->GetNumPoints();
2456  int nqtot = nquad0*nquad1*nquad2;
2457 
2473 
2474  StdRegions::VarCoeffMap::const_iterator MFdir;
2475 
2476  Array<OneD, NekDouble> nFacecdotMF(nqtot,0.0);
2477  Array<OneD, NekDouble> tmp(nqtot);
2478  Array<OneD, NekDouble> tmp_f(nquad_f);
2479  for (int k=0; k<coordim; k++)
2480  {
2481  MFdir = varcoeffs.find(MMFCoeffs[dir*5+k]);
2482  tmp = MFdir->second;
2483 
2484  GetPhysFaceVarCoeffsFromElement(face, FaceExp_f, tmp, tmp_f);
2485 
2486  Vmath::Vvtvp(nquad_f, &tmp_f[0], 1,
2487  &normals[k][0], 1,
2488  &nFacecdotMF[0], 1,
2489  &nFacecdotMF[0], 1);
2490  }
2491 
2492  return nFacecdotMF;
2493  }
void GetPhysFaceVarCoeffsFromElement(const int face, ExpansionSharedPtr &FaceExp, const Array< OneD, const NekDouble > &varcoeff, Array< OneD, NekDouble > &outarray)

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, and Vmath::Vvtvp().

◆ v_GetTraceNormal()

const NormalVector & Nektar::LocalRegions::Expansion3D::v_GetTraceNormal ( const int  face) const
protectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 2495 of file Expansion3D.cpp.

2496  {
2497  auto x = m_faceNormals.find(face);
2498  ASSERTL0 (x != m_faceNormals.end(),
2499  "face normal not computed.");
2500  return x->second;
2501  }
std::map< int, NormalVector > m_faceNormals
Definition: Expansion3D.h:123

References ASSERTL0.

◆ v_GetTraceOrient()

StdRegions::Orientation Nektar::LocalRegions::Expansion3D::v_GetTraceOrient ( int  face)
protectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 2230 of file Expansion3D.cpp.

2231  {
2232  return m_geom->GetForient(face);
2233  }

◆ v_GetTracePhysVals()

void Nektar::LocalRegions::Expansion3D::v_GetTracePhysVals ( const int  face,
const StdRegions::StdExpansionSharedPtr FaceExp,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
StdRegions::Orientation  orient 
)
protectedvirtual

Extract the physical values along face face from inarray into outarray following the local face orientation and point distribution defined by defined in FaceExp.

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 2240 of file Expansion3D.cpp.

2246  {
2247  if (orient == StdRegions::eNoOrientation)
2248  {
2249  orient = GetTraceOrient(face);
2250  }
2251 
2252  int nq0 = FaceExp->GetNumPoints(0);
2253  int nq1 = FaceExp->GetNumPoints(1);
2254 
2255  int nfacepts = GetTraceNumPoints(face);
2256  int dir0 = GetGeom3D()->GetDir(face,0);
2257  int dir1 = GetGeom3D()->GetDir(face,1);
2258 
2259  Array<OneD,NekDouble> o_tmp (nfacepts);
2260  Array<OneD,NekDouble> o_tmp2(FaceExp->GetTotPoints());
2261  Array<OneD, int> faceids;
2262 
2263  // Get local face pts and put into o_tmp
2264  GetTracePhysMap(face,faceids);
2265  // The static cast is necessary because faceids should be
2266  // Array<OneD, size_t> faceids ... or at leaste the same type as
2267  // faceids.size() ...
2268  Vmath::Gathr(static_cast<int>(faceids.size()),inarray,faceids,o_tmp);
2269 
2270  int to_id0,to_id1;
2271 
2273  {
2274  to_id0 = 0;
2275  to_id1 = 1;
2276  }
2277  else // transpose points key evaluation
2278  {
2279  to_id0 = 1;
2280  to_id1 = 0;
2281  }
2282 
2283  // interpolate to points distrbution given in FaceExp
2284  LibUtilities::Interp2D(m_base[dir0]->GetPointsKey(),
2285  m_base[dir1]->GetPointsKey(),
2286  o_tmp.get(),
2287  FaceExp->GetBasis(to_id0)->GetPointsKey(),
2288  FaceExp->GetBasis(to_id1)->GetPointsKey(),
2289  o_tmp2.get());
2290 
2291  // Reshuffule points as required and put into outarray.
2292  v_ReOrientTracePhysMap(orient,faceids, nq0,nq1);
2293  Vmath::Scatr(nq0*nq1,o_tmp2,faceids,outarray);
2294  }
void v_ReOrientTracePhysMap(const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1)
void GetTracePhysMap(const int edge, Array< OneD, int > &outarray)
Definition: Expansion.h:206
int GetTraceNumPoints(const int i) const
This function returns the number of quadrature points belonging to the i-th trace.
Definition: StdExpansion.h:286
void Interp2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
this function interpolates a 2D function evaluated at the quadrature points of the 2D basis,...
Definition: Interp.cpp:115
void Scatr(int n, const T *x, const int *y, T *z)
Scatter vector z[y[i]] = x[i].
Definition: Vmath.cpp:772
void Gathr(int n, const T *sign, const T *x, const int *y, T *z)
Gather vector z[i] = sign[i]*x[y[i]].
Definition: Vmath.cpp:756

References Nektar::StdRegions::eDir1FwdDir2_Dir2FwdDir1, Nektar::StdRegions::eNoOrientation, Vmath::Gathr(), Nektar::LibUtilities::Interp2D(), and Vmath::Scatr().

◆ v_NormVectorIProductWRTBase()

void Nektar::LocalRegions::Expansion3D::v_NormVectorIProductWRTBase ( const Array< OneD, const Array< OneD, NekDouble > > &  Fvec,
Array< OneD, NekDouble > &  outarray 
)
virtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2434 of file Expansion3D.cpp.

2437  {
2438  NormVectorIProductWRTBase(Fvec[0], Fvec[1], Fvec[2], outarray);
2439  }
void NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:628

◆ v_ReOrientTracePhysMap()

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

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 2296 of file Expansion3D.cpp.

2300  {
2301 
2302  if(idmap.size() != nq0*nq1)
2303  {
2304  idmap = Array<OneD,int>(nq0*nq1);
2305  }
2306 
2307  if(GetNverts() == 3) // Tri face
2308  {
2309  switch(orient)
2310  {
2312  // eseentially straight copy
2313  for(int i = 0; i < nq0*nq1; ++i)
2314  {
2315  idmap[i] = i;
2316  }
2317  break;
2319  // reverse
2320  for (int j = 0; j < nq1; ++j)
2321  {
2322  for(int i = 0; i < nq0; ++i)
2323  {
2324  idmap[j*nq0+i] = nq0-1-i +j*nq0;
2325  }
2326  }
2327  break;
2328  default:
2329  ASSERTL0(false,"Case not supposed to be used in this function");
2330  }
2331  }
2332  else
2333  {
2334  switch(orient)
2335  {
2337  // eseentially straight copy
2338  for(int i = 0; i < nq0*nq1; ++i)
2339  {
2340  idmap[i] = i;
2341  }
2342  break;
2344  {
2345  //Direction A negative and B positive
2346  for (int j = 0; j < nq1; j++)
2347  {
2348  for (int i =0; i < nq0; ++i)
2349  {
2350  idmap[j*nq0+i] = nq0-1-i + j*nq0;
2351  }
2352  }
2353  }
2354  break;
2356  {
2357  //Direction A positive and B negative
2358  for (int j = 0; j < nq1; j++)
2359  {
2360  for (int i =0; i < nq0; ++i)
2361  {
2362  idmap[j*nq0+i] = nq0*(nq1-1-j)+i;
2363  }
2364  }
2365  }
2366  break;
2368  {
2369  //Direction A negative and B negative
2370  for (int j = 0; j < nq1; j++)
2371  {
2372  for (int i =0; i < nq0; ++i)
2373  {
2374  idmap[j*nq0+i] = nq0*nq1-1-j*nq0 -i;
2375  }
2376  }
2377  }
2378  break;
2380  {
2381  //Transposed, Direction A and B positive
2382  for (int i =0; i < nq0; ++i)
2383  {
2384  for (int j = 0; j < nq1; ++j)
2385  {
2386  idmap[i*nq1+j] = i + j*nq0;
2387  }
2388  }
2389  }
2390  break;
2392  {
2393  //Transposed, Direction A positive and B negative
2394  for (int i =0; i < nq0; ++i)
2395  {
2396  for (int j = 0; j < nq1; ++j)
2397  {
2398  idmap[i*nq1+j] = nq0-1-i + j*nq0;
2399  }
2400  }
2401  }
2402  break;
2404  {
2405  //Transposed, Direction A negative and B positive
2406  for (int i =0; i < nq0; ++i)
2407  {
2408  for (int j = 0; j < nq1; ++j)
2409  {
2410  idmap[i*nq1+j] = i+nq0*(nq1-1)-j*nq0;
2411  }
2412  }
2413  }
2414  break;
2416  {
2417  //Transposed, Direction A and B negative
2418  for (int i =0; i < nq0; ++i)
2419  {
2420  for (int j = 0; j < nq1; ++j)
2421  {
2422  idmap[i*nq1+j] = nq0*nq1-1-i-j*nq0;
2423  }
2424  }
2425  }
2426  break;
2427  default:
2428  ASSERTL0(false,"Unknow orientation");
2429  break;
2430  }
2431  }
2432  }

References ASSERTL0, Nektar::StdRegions::eDir1BwdDir1_Dir2BwdDir2, Nektar::StdRegions::eDir1BwdDir1_Dir2FwdDir2, Nektar::StdRegions::eDir1BwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1BwdDir2_Dir2FwdDir1, Nektar::StdRegions::eDir1FwdDir1_Dir2BwdDir2, Nektar::StdRegions::eDir1FwdDir1_Dir2FwdDir2, Nektar::StdRegions::eDir1FwdDir2_Dir2BwdDir1, and Nektar::StdRegions::eDir1FwdDir2_Dir2FwdDir1.

Member Data Documentation

◆ m_faceNormals

std::map<int, NormalVector> Nektar::LocalRegions::Expansion3D::m_faceNormals
protected

◆ m_requireNeg

std::vector<bool> Nektar::LocalRegions::Expansion3D::m_requireNeg
private

Definition at line 188 of file Expansion3D.h.