Nektar++
Loading...
Searching...
No Matches
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
Nektar::LocalRegions::Expansion Class Reference

#include <Expansion.h>

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

Public Member Functions

 Expansion (SpatialDomains::Geometry *pGeom)
 
 Expansion (const Expansion &pSrc)
 
 ~Expansion () override
 
void SetTraceExp (const int traceid, ExpansionSharedPtr &f)
 
ExpansionSharedPtr GetTraceExp (const int traceid)
 
DNekScalMatSharedPtr GetLocMatrix (const LocalRegions::MatrixKey &mkey)
 
void DropLocMatrix (const LocalRegions::MatrixKey &mkey)
 
DNekScalMatSharedPtr GetLocMatrix (const StdRegions::MatrixType mtype, const StdRegions::ConstFactorMap &factors=StdRegions::NullConstFactorMap, const StdRegions::VarCoeffMap &varcoeffs=StdRegions::NullVarCoeffMap)
 
SpatialDomains::GeometryGetGeom () const
 
void Reset ()
 
IndexMapValuesSharedPtr CreateIndexMap (const IndexMapKey &ikey)
 
DNekScalBlkMatSharedPtr CreateStaticCondMatrix (const MatrixKey &mkey)
 
const SpatialDomains::GeomFactorsSharedPtrGetMetricInfo () const
 
DNekMatSharedPtr BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
 
DNekMatSharedPtr BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd)
 
void ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType)
 
void AddEdgeNormBoundaryInt (const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
void AddEdgeNormBoundaryInt (const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
void AddFaceNormBoundaryInt (const int face, const std::shared_ptr< Expansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
void DGDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &coeffs, Array< OneD, NekDouble > &outarray)
 
NekDouble VectorFlux (const Array< OneD, Array< OneD, NekDouble > > &vec)
 
void NormalTraceDerivFactors (Array< OneD, Array< OneD, NekDouble > > &factors, Array< OneD, Array< OneD, NekDouble > > &d0factors, Array< OneD, Array< OneD, NekDouble > > &d1factors)
 
IndexMapValuesSharedPtr GetIndexMap (const IndexMapKey &ikey)
 
void AlignVectorToCollapsedDir (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
ExpansionSharedPtr GetLeftAdjacentElementExp () const
 
ExpansionSharedPtr GetRightAdjacentElementExp () const
 
int GetLeftAdjacentElementTrace () const
 
int GetRightAdjacentElementTrace () const
 
void SetAdjacentElementExp (int traceid, ExpansionSharedPtr &e)
 
StdRegions::Orientation GetTraceOrient (int trace)
 
void SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void DivideByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Divided by the metric jacobi and quadrature weights.
 
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).
 
void GetTracePhysVals (const int trace, const StdRegions::StdExpansionSharedPtr &TraceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=StdRegions::eNoOrientation)
 
void GetTracePhysMap (const int edge, Array< OneD, int > &outarray)
 
void ReOrientTracePhysMap (const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1)
 
const NormalVectorGetTraceNormal (const int id)
 
void ComputeTraceNormal (const int id)
 
const Array< OneD, const NekDouble > & GetPhysNormals (void)
 
void SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
void SetUpPhysNormals (const int trace)
 
void AddRobinMassMatrix (const int traceid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
void TraceNormLen (const int traceid, NekDouble &h, NekDouble &p)
 
void AddRobinTraceContribution (const int traceid, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs)
 
const Array< OneD, const NekDouble > & GetElmtBndNormDirElmtLen (const int nbnd) const
 
void StdDerivBaseOnTraceMat (Array< OneD, DNekMatSharedPtr > &DerivMat)
 
void PhysDerivBaseOnTraceMat (const int traceid, Array< OneD, DNekMatSharedPtr > &DerivMat)
 
void PhysBaseOnTraceMat (const int traceid, DNekMatSharedPtr &BdataMat)
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion
 StdExpansion ()
 Default Constructor.
 
 StdExpansion (const int numcoeffs, const int numbases, const LibUtilities::BasisKey &Ba=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bb=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bc=LibUtilities::NullBasisKey)
 Constructor.
 
 StdExpansion (const StdExpansion &T)
 Copy Constructor.
 
virtual ~StdExpansion ()
 Destructor.
 
int GetNumBases () const
 This function returns the number of 1D bases used in the expansion.
 
const Array< OneD, const LibUtilities::BasisSharedPtr > & GetBase () const
 This function gets the shared point to basis.
 
const LibUtilities::BasisSharedPtrGetBasis (int dir) const
 This function gets the shared point to basis in the dir direction.
 
int GetNcoeffs (void) const
 This function returns the total number of coefficients used in the expansion.
 
int GetTotPoints () const
 This function returns the total number of quadrature points used in the element.
 
LibUtilities::BasisType GetBasisType (const int dir) const
 This function returns the type of basis used in the dir direction.
 
int GetBasisNumModes (const int dir) const
 This function returns the number of expansion modes in the dir direction.
 
int EvalBasisNumModesMax (void) const
 This function returns the maximum number of expansion modes over all local directions.
 
LibUtilities::PointsType GetPointsType (const int dir) const
 This function returns the type of quadrature points used in the dir direction.
 
int GetNumPoints (const int dir) const
 This function returns the number of quadrature points in the dir direction.
 
const Array< OneD, const NekDouble > & GetPoints (const int dir) const
 This function returns a pointer to the array containing the quadrature points in dir direction.
 
int GetNverts () const
 This function returns the number of vertices of the expansion domain.
 
int GetTraceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th trace.
 
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.
 
const LibUtilities::BasisKey GetTraceBasisKey (const int i, int k=-1, bool UseGLL=false) const
 This function returns the basis key belonging to the i-th trace.
 
LibUtilities::PointsKey GetTracePointsKey (const int i, int k=-1) const
 This function returns the basis key belonging to the i-th trace.
 
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.
 
int GetNtraces () const
 Returns the number of trace elements connected to this element.
 
LibUtilities::ShapeType DetShapeType () const
 This function returns the shape of the expansion domain.
 
std::shared_ptr< StdExpansionGetStdExp () const
 
std::shared_ptr< StdExpansionGetLinStdExp (void) const
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion () const
 
bool IsNodalNonTensorialExp ()
 
void NodalToModal (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space.
 
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space.
 
void FwdTransBndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble Integral (const Array< OneD, const NekDouble > &inarray)
 This function integrates the specified function over the domain.
 
void FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 This function fills the array outarray with the mode-th mode of the expansion.
 
void IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 this function calculates the inner product of a given function f with the different modes of the expansion
 
void IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check)
 
void IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
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.
 
void SetElmtId (const int id)
 Set the element id of this expansion when used in a list by returning value of m_elmt_id.
 
void GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2=NullNekDouble1DArray, Array< OneD, NekDouble > &coords_3=NullNekDouble1DArray)
 this function returns the physical coordinates of the quadrature points of the expansion
 
void GetCoord (const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
 given the coordinates of a point of the element in the local collapsed coordinate system, this function calculates the physical coordinates of the point
 
DNekMatSharedPtr GetStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr GetStdStaticCondMatrix (const StdMatrixKey &mkey)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
void DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
int CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
NekDouble StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
int GetCoordim ()
 
void GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
void GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
int GetVertexMap (const int localVertexId, bool useCoeffPacking=false)
 
void GetTraceToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
 
void GetTraceCoeffMap (const unsigned int traceid, Array< OneD, unsigned int > &maparray)
 
void GetElmtTraceToTraceMap (const unsigned int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
 
void GetTraceInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eForwards)
 
void GetTraceNumModes (const int tid, int &numModes0, int &numModes1, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2)
 
void MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
DNekMatSharedPtr CreateGeneralMatrix (const StdMatrixKey &mkey)
 this function generates the mass matrix \(\mathbf{M}[i][j] = \int \phi_i(\mathbf{x}) \phi_j(\mathbf{x}) d\mathbf{x}\)
 
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 LinearAdvectionMatrixOp (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)
 
NekDouble PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain.
 
NekDouble PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
 This function evaluates the first derivative of the expansion at a single (arbitrary) point of the domain.
 
NekDouble PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs, std::array< NekDouble, 6 > &secondOrderDerivs)
 
NekDouble PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain.
 
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.
 
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta.
 
void LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi)
 Convert local collapsed coordinates eta into local cartesian coordinate xi.
 
void PhysInterp (std::shared_ptr< StdExpansion > fromExp, const Array< OneD, const NekDouble > &fromData, Array< OneD, NekDouble > &toData)
 interpolate from one set of quadrature points available from FromExp to the set of quadrature points in the current expansion. If the points are the same this routine will just copy the data
 
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual void v_DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
NekDouble Linf (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( L_\infty\) error \( |\epsilon|_\infty = \max |u - u_{exact}|\) where \( u_{exact}\) is given by the array sol.
 
NekDouble L2 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( L_2\) error, \( | \epsilon |_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 dx \right]^{1/2} d\xi_1 \) where \( u_{exact}\) is given by the array sol.
 
NekDouble H1 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( H^1\) error, \( | \epsilon |^1_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 + \nabla(u - u_{exact})\cdot\nabla(u - u_{exact})\cdot dx \right]^{1/2} d\xi_1 \) where \( u_{exact}\) is given by the array sol.
 
const 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.
 
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.
 
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.
 
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)
 

Protected Member Functions

void ComputeLaplacianMetric ()
 
void ComputeQuadratureMetric ()
 
void ComputeGmatcdotMF (const Array< TwoD, const NekDouble > &df, const Array< OneD, const NekDouble > &direction, Array< OneD, Array< OneD, NekDouble > > &dfdir)
 
Array< OneD, NekDoubleGetMF (const int dir, const int shapedim, const StdRegions::VarCoeffMap &varcoeffs)
 
Array< OneD, NekDoubleGetMFDiv (const int dir, const StdRegions::VarCoeffMap &varcoeffs)
 
Array< OneD, NekDoubleGetMFMag (const int dir, const StdRegions::VarCoeffMap &varcoeffs)
 
void v_MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
virtual void v_DivideByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_ComputeLaplacianMetric ()
 
int v_GetCoordim () const override
 
void v_GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
 
virtual DNekScalMatSharedPtr v_GetLocMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual void v_DropLocMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual DNekMatSharedPtr v_BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
 
virtual DNekMatSharedPtr v_BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd)
 
virtual void v_ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType)
 
virtual void v_AddEdgeNormBoundaryInt (const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddEdgeNormBoundaryInt (const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddFaceNormBoundaryInt (const int face, const std::shared_ptr< Expansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
virtual void v_DGDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &coeffs, Array< OneD, NekDouble > &outarray)
 
virtual NekDouble v_VectorFlux (const Array< OneD, Array< OneD, NekDouble > > &vec)
 
virtual void v_NormalTraceDerivFactors (Array< OneD, Array< OneD, NekDouble > > &factors, Array< OneD, Array< OneD, NekDouble > > &d0factors, Array< OneD, Array< OneD, NekDouble > > &d1factors)
 
virtual void v_AlignVectorToCollapsedDir (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual StdRegions::Orientation v_GetTraceOrient (int trace)
 
void v_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
virtual void v_GetTraceQFactors (const int trace, Array< OneD, NekDouble > &outarray)
 
virtual void v_GetTracePhysVals (const int trace, const StdRegions::StdExpansionSharedPtr &TraceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient)
 
virtual void v_GetTracePhysMap (const int edge, Array< OneD, int > &outarray)
 
virtual void v_ReOrientTracePhysMap (const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1=-1)
 
virtual void v_ComputeTraceNormal (const int id)
 
virtual const Array< OneD, const NekDouble > & v_GetPhysNormals ()
 
virtual void v_SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
virtual void v_SetUpPhysNormals (const int id)
 
virtual void v_AddRobinMassMatrix (const int face, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
virtual void v_AddRobinTraceContribution (const int traceid, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs)
 
virtual void v_TraceNormLen (const int traceid, NekDouble &h, NekDouble &p)
 
virtual void v_GenTraceExp (const int traceid, ExpansionSharedPtr &exp)
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion
DNekMatSharedPtr CreateStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr CreateStdStaticCondMatrix (const StdMatrixKey &mkey)
 Create the static condensation of a matrix when using a boundary interior decomposition.
 
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 LinearAdvectionMatrixOp_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 NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
virtual void v_GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 
virtual void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals, NekDouble &deriv, NekDouble &deriv2)
 This function performs the barycentric interpolation of the polynomial stored in coord at a point physvals using barycentric interpolation weights in direction.
 
template<int DIR>
NekDouble BaryEvaluateBasis (const NekDouble &coord, const int &mode)
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals)
 Helper function to pass an unused value by reference into BaryEvaluate.
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals, NekDouble &deriv)
 

Protected Attributes

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

Detailed Description

Definition at line 71 of file Expansion.h.

Constructor & Destructor Documentation

◆ Expansion() [1/2]

Nektar::LocalRegions::Expansion::Expansion ( SpatialDomains::Geometry pGeom)

Definition at line 43 of file Expansion.cpp.

45 std::bind(&Expansion::CreateIndexMap, this, std::placeholders::_1),
46 std::string("ExpansionIndexMap")),
49{
50 if (!m_metricinfo)
51 {
52 return;
53 }
54
55 if (!m_metricinfo->IsValid())
56 {
57 int nDim = m_base.size();
58 string type = "regular";
59 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
60 {
61 type = "deformed";
62 }
63
64 stringstream err;
65 err << nDim << "D " << type << " Jacobian not positive "
66 << "(element ID = " << m_geom->GetGlobalID() << ") "
67 << "(first vertex ID = " << m_geom->GetVid(0) << ")";
68 NEKERROR(ErrorUtil::ewarning, err.str());
69 }
70
71 m_traceExp.clear();
72}
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
IndexMapValuesSharedPtr CreateIndexMap(const IndexMapKey &ikey)
SpatialDomains::Geometry * m_geom
Definition Expansion.h:279
std::map< int, ExpansionWeakPtr > m_traceExp
Definition Expansion.h:278
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition Expansion.h:280
LibUtilities::NekManager< IndexMapKey, IndexMapValues, IndexMapKey::opLess > m_indexMapManager
Definition Expansion.h:276
int GetVid(int i) const
Returns global id of vertex i of this object.
Definition Geometry.h:353
int GetGlobalID(void) const
Get the ID of this object.
Definition Geometry.h:322
GeomFactorsSharedPtr GetGeomFactors()
Get the geometric factors for this object, generating them if required.
Definition Geometry.h:297
Array< OneD, LibUtilities::BasisSharedPtr > m_base
@ eDeformed
Geometry is curved or has non-constant factors.

References Nektar::SpatialDomains::eDeformed, Nektar::ErrorUtil::ewarning, Nektar::SpatialDomains::Geometry::GetGlobalID(), Nektar::SpatialDomains::Geometry::GetVid(), Nektar::StdRegions::StdExpansion::m_base, m_geom, m_metricinfo, m_traceExp, and NEKERROR.

◆ Expansion() [2/2]

Nektar::LocalRegions::Expansion::Expansion ( const Expansion pSrc)

Definition at line 74 of file Expansion.cpp.

75 : StdExpansion(pSrc), m_indexMapManager(pSrc.m_indexMapManager),
76 m_geom(pSrc.m_geom), m_metricinfo(pSrc.m_metricinfo)
77{
78}
StdExpansion()
Default Constructor.

◆ ~Expansion()

Nektar::LocalRegions::Expansion::~Expansion ( )
override

Definition at line 80 of file Expansion.cpp.

81{
82}

Member Function Documentation

◆ AddEdgeNormBoundaryInt() [1/2]

void Nektar::LocalRegions::Expansion::AddEdgeNormBoundaryInt ( const int  edge,
const std::shared_ptr< Expansion > &  EdgeExp,
const Array< OneD, const NekDouble > &  Fn,
Array< OneD, NekDouble > &  outarray 
)

Definition at line 122 of file Expansion.cpp.

125{
126 v_AddEdgeNormBoundaryInt(edge, EdgeExp, Fn, outarray);
127}
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)

References v_AddEdgeNormBoundaryInt().

◆ AddEdgeNormBoundaryInt() [2/2]

void Nektar::LocalRegions::Expansion::AddEdgeNormBoundaryInt ( const int  edge,
const std::shared_ptr< Expansion > &  EdgeExp,
const Array< OneD, const NekDouble > &  Fx,
const Array< OneD, const NekDouble > &  Fy,
Array< OneD, NekDouble > &  outarray 
)

Definition at line 114 of file Expansion.cpp.

118{
119 v_AddEdgeNormBoundaryInt(edge, EdgeExp, Fx, Fy, outarray);
120}

References v_AddEdgeNormBoundaryInt().

Referenced by Nektar::LocalRegions::Expansion2D::v_AddEdgeNormBoundaryInt().

◆ AddFaceNormBoundaryInt()

void Nektar::LocalRegions::Expansion::AddFaceNormBoundaryInt ( const int  face,
const std::shared_ptr< Expansion > &  FaceExp,
const Array< OneD, const NekDouble > &  Fn,
Array< OneD, NekDouble > &  outarray 
)

Definition at line 129 of file Expansion.cpp.

132{
133 v_AddFaceNormBoundaryInt(face, FaceExp, Fn, outarray);
134}
virtual void v_AddFaceNormBoundaryInt(const int face, const std::shared_ptr< Expansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)

References v_AddFaceNormBoundaryInt().

◆ AddRobinMassMatrix()

void Nektar::LocalRegions::Expansion::AddRobinMassMatrix ( const int  traceid,
const Array< OneD, const NekDouble > &  primCoeffs,
DNekMatSharedPtr inoutmat 
)
inline

Definition at line 243 of file Expansion.h.

246 {
247 v_AddRobinMassMatrix(traceid, primCoeffs, inoutmat);
248 }
virtual void v_AddRobinMassMatrix(const int face, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)

References v_AddRobinMassMatrix().

◆ AddRobinTraceContribution()

void Nektar::LocalRegions::Expansion::AddRobinTraceContribution ( const int  traceid,
const Array< OneD, const NekDouble > &  primCoeffs,
const Array< OneD, NekDouble > &  incoeffs,
Array< OneD, NekDouble > &  coeffs 
)
inline

Definition at line 255 of file Expansion.h.

258 {
259 v_AddRobinTraceContribution(traceid, primCoeffs, incoeffs, coeffs);
260 }
virtual void v_AddRobinTraceContribution(const int traceid, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs)

References v_AddRobinTraceContribution().

◆ AlignVectorToCollapsedDir()

void Nektar::LocalRegions::Expansion::AlignVectorToCollapsedDir ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
inline

Definition at line 151 of file Expansion.h.

154 {
155 v_AlignVectorToCollapsedDir(dir, inarray, outarray);
156 }
virtual void v_AlignVectorToCollapsedDir(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)

References v_AlignVectorToCollapsedDir().

Referenced by Nektar::LocalRegions::NodalTriExp::v_IProductWRTDerivBase_SumFac().

◆ BuildTransformationMatrix()

DNekMatSharedPtr Nektar::LocalRegions::Expansion::BuildTransformationMatrix ( const DNekScalMatSharedPtr r_bnd,
const StdRegions::MatrixType  matrixType 
)

Definition at line 95 of file Expansion.cpp.

97{
98 return v_BuildTransformationMatrix(r_bnd, matrixType);
99}
virtual DNekMatSharedPtr v_BuildTransformationMatrix(const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)

References v_BuildTransformationMatrix().

Referenced by Nektar::LocalRegions::Expansion3D::CreateMatrix().

◆ BuildVertexMatrix()

DNekMatSharedPtr Nektar::LocalRegions::Expansion::BuildVertexMatrix ( const DNekScalMatSharedPtr r_bnd)

Definition at line 101 of file Expansion.cpp.

102{
103 return v_BuildVertexMatrix(r_bnd);
104}
virtual DNekMatSharedPtr v_BuildVertexMatrix(const DNekScalMatSharedPtr &r_bnd)

References v_BuildVertexMatrix().

Referenced by Nektar::LocalRegions::Expansion2D::CreateMatrix(), and Nektar::LocalRegions::Expansion3D::CreateMatrix().

◆ ComputeGmatcdotMF()

void Nektar::LocalRegions::Expansion::ComputeGmatcdotMF ( const Array< TwoD, const NekDouble > &  df,
const Array< OneD, const NekDouble > &  direction,
Array< OneD, Array< OneD, NekDouble > > &  dfdir 
)
protected

Definition at line 688 of file Expansion.cpp.

691{
692 int shapedim = dfdir.size();
693 int coordim = GetCoordim();
694 int nqtot = direction.size() / coordim;
695
696 for (int j = 0; j < shapedim; j++)
697 {
698 dfdir[j] = Array<OneD, NekDouble>(nqtot, 0.0);
699 for (int k = 0; k < coordim; k++)
700 {
701 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
702 {
703 Vmath::Vvtvp(nqtot, &df[shapedim * k + j][0], 1,
704 &direction[k * nqtot], 1, &dfdir[j][0], 1,
705 &dfdir[j][0], 1);
706 }
707 else
708 {
709 Vmath::Svtvp(nqtot, df[shapedim * k + j][0],
710 &direction[k * nqtot], 1, &dfdir[j][0], 1,
711 &dfdir[j][0], 1);
712 }
713 }
714 }
715}
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.hpp:396
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.hpp:366

References Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetCoordim(), m_metricinfo, Vmath::Svtvp(), and Vmath::Vvtvp().

Referenced by Nektar::LocalRegions::Expansion2D::CreateMatrix(), Nektar::LocalRegions::HexExp::v_IProductWRTDirectionalDerivBase_SumFac(), Nektar::LocalRegions::TriExp::v_IProductWRTDirectionalDerivBase_SumFac(), and Nektar::LocalRegions::HexExp::v_PhysDirectionalDeriv().

◆ ComputeLaplacianMetric()

void Nektar::LocalRegions::Expansion::ComputeLaplacianMetric ( )
protected

◆ ComputeQuadratureMetric()

void Nektar::LocalRegions::Expansion::ComputeQuadratureMetric ( )
protected

Definition at line 466 of file Expansion.cpp.

467{
468 unsigned int nqtot = GetTotPoints();
469 SpatialDomains::GeomType type = m_metricinfo->GetGtype();
471 if (type == SpatialDomains::eRegular ||
473 {
475 Array<OneD, NekDouble>(nqtot, m_metricinfo->GetJac(p)[0]);
476 }
477 else
478 {
480 }
481
484}
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
const LibUtilities::PointsKeyVector GetPointsKeys() const
virtual void v_MultiplyByStdQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
std::vector< PointsKey > PointsKeyVector
Definition Points.h:231
GeomType
Indicates the type of element geometry.
@ eRegular
Geometry is straight-sided with constant geometric factors.
@ eMovingRegular
Currently unused.
std::vector< double > p(NPUPPER)

References Nektar::LocalRegions::eMetricQuadrature, Nektar::SpatialDomains::eMovingRegular, Nektar::SpatialDomains::eRegular, Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::GetTotPoints(), m_metricinfo, m_metrics, and Nektar::StdRegions::StdExpansion::v_MultiplyByStdQuadratureMetric().

Referenced by Nektar::LocalRegions::HexExp::v_ComputeLaplacianMetric(), Nektar::LocalRegions::PyrExp::v_ComputeLaplacianMetric(), Nektar::LocalRegions::QuadExp::v_ComputeLaplacianMetric(), Nektar::LocalRegions::TriExp::v_ComputeLaplacianMetric(), v_DivideByQuadratureMetric(), and v_MultiplyByQuadratureMetric().

◆ ComputeTraceNormal()

void Nektar::LocalRegions::Expansion::ComputeTraceNormal ( const int  id)
inline

Definition at line 223 of file Expansion.h.

224 {
226 }
virtual void v_ComputeTraceNormal(const int id)

References v_ComputeTraceNormal().

◆ CreateIndexMap()

IndexMapValuesSharedPtr Nektar::LocalRegions::Expansion::CreateIndexMap ( const IndexMapKey ikey)

Definition at line 181 of file Expansion.cpp.

182{
183 IndexMapValuesSharedPtr returnval;
184
185 IndexMapType itype = ikey.GetIndexMapType();
186
187 int entity = ikey.GetIndexEntity();
188
189 StdRegions::Orientation orient = ikey.GetIndexOrientation();
190
191 Array<OneD, unsigned int> map;
192 Array<OneD, int> sign;
193
194 switch (itype)
195 {
196 case eEdgeToElement:
197 {
198 GetTraceToElementMap(entity, map, sign, orient);
199 }
200 break;
201 case eFaceToElement:
202 {
203 GetTraceToElementMap(entity, map, sign, orient);
204 }
205 break;
206 case eEdgeInterior:
207 {
208 ASSERTL0(false, "Boundary Index Map not implemented yet.");
209 // v_GetEdgeInteriorMap(entity,orient,map,sign);
210 }
211 break;
212 case eFaceInterior:
213 {
214 ASSERTL0(false, "Boundary Index Map not implemented yet.");
215 // v_GetFaceInteriorMap(entity,orient,map,sign);
216 }
217 break;
218 case eBoundary:
219 {
220 ASSERTL0(false, "Boundary Index Map not implemented yet.");
221 }
222 break;
223 case eVertex:
224 {
225 ASSERTL0(false, "Vertex Index Map not implemented yet.");
226 }
227 break;
228 default:
229 {
230 ASSERTL0(false, "The Index Map you are requiring "
231 "is not between the possible options.");
232 }
233 }
234
236
237 for (int i = 0; i < map.size(); i++)
238 {
239 (*returnval)[i].index = map[i];
240 (*returnval)[i].sign = sign[i];
241 }
242
243 return returnval;
244}
#define ASSERTL0(condition, msg)
#define sign(a, b)
return the sign(b)*a
Definition Polylib.cpp:47
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void GetTraceToElementMap(const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
std::shared_ptr< IndexMapValues > IndexMapValuesSharedPtr

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::LocalRegions::eBoundary, Nektar::LocalRegions::eEdgeInterior, Nektar::LocalRegions::eEdgeToElement, Nektar::LocalRegions::eFaceInterior, Nektar::LocalRegions::eFaceToElement, Nektar::LocalRegions::eVertex, Nektar::LocalRegions::IndexMapKey::GetIndexEntity(), Nektar::LocalRegions::IndexMapKey::GetIndexMapType(), Nektar::LocalRegions::IndexMapKey::GetIndexOrientation(), Nektar::StdRegions::StdExpansion::GetTraceToElementMap(), and sign.

◆ CreateStaticCondMatrix()

DNekScalBlkMatSharedPtr Nektar::LocalRegions::Expansion::CreateStaticCondMatrix ( const MatrixKey mkey)

Definition at line 272 of file Expansion.cpp.

273{
274 DNekScalBlkMatSharedPtr returnval;
275
277 "Geometric information is not set up");
278
279 // set up block matrix system
280 unsigned int nbdry = NumBndryCoeffs();
281 unsigned int nint = (unsigned int)(m_ncoeffs - nbdry);
282 unsigned int exp_size[] = {nbdry, nint};
283 unsigned int nblks = 2;
285 nblks, nblks, exp_size, exp_size);
286 // Really need a constructor which takes Arrays
287 NekDouble factor = 1.0;
288
289 switch (mkey.GetMatrixType())
290 {
291 // this can only use stdregions statically condensed system
292 // for mass matrix
294 if ((m_metricinfo->GetGtype() == SpatialDomains::eDeformed) ||
295 (mkey.GetNVarCoeff()))
296 {
297 factor = 1.0;
298 goto UseLocRegionsMatrix;
299 }
300 else
301 {
302 factor = (m_metricinfo->GetJac(GetPointsKeys()))[0];
303 goto UseStdRegionsMatrix;
304 }
305 break;
306 default: // use Deformed case for both
307 // regular and deformed geometries
308 factor = 1.0;
309 goto UseLocRegionsMatrix;
310 break;
311 UseStdRegionsMatrix:
312 {
313 NekDouble invfactor = 1.0 / factor;
314 NekDouble one = 1.0;
317 DNekMatSharedPtr Asubmat;
318
319 returnval->SetBlock(
320 0, 0,
322 factor, Asubmat = mat->GetBlock(0, 0)));
323 returnval->SetBlock(
324 0, 1,
326 one, Asubmat = mat->GetBlock(0, 1)));
327 returnval->SetBlock(
328 1, 0,
330 factor, Asubmat = mat->GetBlock(1, 0)));
331 returnval->SetBlock(
332 1, 1,
334 invfactor, Asubmat = mat->GetBlock(1, 1)));
335 }
336 break;
337 UseLocRegionsMatrix:
338 {
339 int i, j;
340 NekDouble invfactor = 1.0 / factor;
341 NekDouble one = 1.0;
342 DNekScalMat &mat = *GetLocMatrix(mkey);
351
352 Array<OneD, unsigned int> bmap(nbdry);
353 Array<OneD, unsigned int> imap(nint);
354 GetBoundaryMap(bmap);
355 GetInteriorMap(imap);
356
357 for (i = 0; i < nbdry; ++i)
358 {
359 for (j = 0; j < nbdry; ++j)
360 {
361 (*A)(i, j) = mat(bmap[i], bmap[j]);
362 }
363
364 for (j = 0; j < nint; ++j)
365 {
366 (*B)(i, j) = mat(bmap[i], imap[j]);
367 }
368 }
369
370 for (i = 0; i < nint; ++i)
371 {
372 for (j = 0; j < nbdry; ++j)
373 {
374 (*C)(i, j) = mat(imap[i], bmap[j]);
375 }
376
377 for (j = 0; j < nint; ++j)
378 {
379 (*D)(i, j) = mat(imap[i], imap[j]);
380 }
381 }
382
383 // Calculate static condensed system
384 if (nint)
385 {
386 D->Invert();
387 (*B) = (*B) * (*D);
388 (*A) = (*A) - (*B) * (*C);
389 }
390
392
393 returnval->SetBlock(
394 0, 0,
395 Atmp =
397 returnval->SetBlock(
398 0, 1,
400 returnval->SetBlock(
401 1, 0,
402 Atmp =
404 returnval->SetBlock(
405 1, 1,
407 D));
408
409 // Remove the local matrix from manager if using this option since
410 // we assume it is only created to generate static condensed system
411 v_DropLocMatrix(mkey);
412 }
413 }
414 return returnval;
415}
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
virtual void v_DropLocMatrix(const LocalRegions::MatrixKey &mkey)
DNekScalMatSharedPtr GetLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition Expansion.cpp:84
void GetBoundaryMap(Array< OneD, unsigned int > &outarray)
DNekBlkMatSharedPtr GetStdStaticCondMatrix(const StdMatrixKey &mkey)
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
@ eNoGeomType
No type defined.
NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag > DNekScalMat
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
std::shared_ptr< DNekMat > DNekMatSharedPtr

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL2, Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::eMass, Nektar::SpatialDomains::eNoGeomType, Nektar::StdRegions::StdExpansion::GetBoundaryMap(), Nektar::StdRegions::StdExpansion::GetInteriorMap(), GetLocMatrix(), Nektar::StdRegions::StdMatrixKey::GetMatrixType(), Nektar::StdRegions::StdMatrixKey::GetNVarCoeff(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::GetStdStaticCondMatrix(), m_metricinfo, Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::StdRegions::StdExpansion::NumBndryCoeffs(), and v_DropLocMatrix().

◆ DGDeriv()

void Nektar::LocalRegions::Expansion::DGDeriv ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, ExpansionSharedPtr > &  EdgeExp,
Array< OneD, Array< OneD, NekDouble > > &  coeffs,
Array< OneD, NekDouble > &  outarray 
)

Definition at line 136 of file Expansion.cpp.

141{
142 v_DGDeriv(dir, inarray, EdgeExp, coeffs, outarray);
143}
virtual void v_DGDeriv(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &coeffs, Array< OneD, NekDouble > &outarray)

References v_DGDeriv().

◆ DivideByQuadratureMetric()

void Nektar::LocalRegions::Expansion::DivideByQuadratureMetric ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
inline

Divided by the metric jacobi and quadrature weights.

Definition at line 181 of file Expansion.h.

184 {
185 v_DivideByQuadratureMetric(inarray, outarray);
186 }
virtual void v_DivideByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)

References v_DivideByQuadratureMetric().

◆ DropLocMatrix()

void Nektar::LocalRegions::Expansion::DropLocMatrix ( const LocalRegions::MatrixKey mkey)

◆ ExtractDataToCoeffs()

void Nektar::LocalRegions::Expansion::ExtractDataToCoeffs ( const NekDouble data,
const std::vector< unsigned int > &  nummodes,
const int  nmodes_offset,
NekDouble coeffs,
std::vector< LibUtilities::BasisType > &  fromType 
)

Definition at line 106 of file Expansion.cpp.

110{
111 v_ExtractDataToCoeffs(data, nummodes, nmodes_offset, coeffs, fromType);
112}
virtual void v_ExtractDataToCoeffs(const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType)

References v_ExtractDataToCoeffs().

Referenced by Nektar::LocalRegions::NodalTriExp::v_ExtractDataToCoeffs().

◆ GetElmtBndNormDirElmtLen()

const Array< OneD, const NekDouble > & Nektar::LocalRegions::Expansion::GetElmtBndNormDirElmtLen ( const int  nbnd) const

Definition at line 999 of file Expansion.cpp.

1001{
1002 auto x = m_elmtBndNormDirElmtLen.find(nbnd);
1004 "m_elmtBndNormDirElmtLen normal not computed.");
1005 return x->second;
1006}
std::map< int, Array< OneD, NekDouble > > m_elmtBndNormDirElmtLen
the element length in each element boundary(Vertex, edge or face) normal direction calculated based o...
Definition Expansion.h:292

References ASSERTL0, and m_elmtBndNormDirElmtLen.

◆ GetGeom()

SpatialDomains::Geometry * Nektar::LocalRegions::Expansion::GetGeom ( ) const

◆ GetIndexMap()

IndexMapValuesSharedPtr Nektar::LocalRegions::Expansion::GetIndexMap ( const IndexMapKey ikey)
inline

◆ GetLeftAdjacentElementExp()

ExpansionSharedPtr Nektar::LocalRegions::Expansion::GetLeftAdjacentElementExp ( ) const
inline

Definition at line 447 of file Expansion.h.

448{
449 ASSERTL1(m_elementLeft.lock().get(), "Left adjacent element not set.");
450 return m_elementLeft.lock();
451}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
ExpansionWeakPtr m_elementLeft
Definition Expansion.h:283

References ASSERTL1, and m_elementLeft.

Referenced by Nektar::LocalRegions::PointExp::v_NormVectorIProductWRTBase(), Nektar::LocalRegions::SegExp::v_NormVectorIProductWRTBase(), Nektar::LocalRegions::QuadExp::v_NormVectorIProductWRTBase(), Nektar::LocalRegions::TriExp::v_NormVectorIProductWRTBase(), Nektar::LocalRegions::Expansion1D::v_VectorFlux(), and Nektar::LocalRegions::Expansion2D::v_VectorFlux().

◆ GetLeftAdjacentElementTrace()

int Nektar::LocalRegions::Expansion::GetLeftAdjacentElementTrace ( ) const
inline

◆ GetLocMatrix() [1/2]

DNekScalMatSharedPtr Nektar::LocalRegions::Expansion::GetLocMatrix ( const LocalRegions::MatrixKey mkey)

◆ GetLocMatrix() [2/2]

DNekScalMatSharedPtr Nektar::LocalRegions::Expansion::GetLocMatrix ( const StdRegions::MatrixType  mtype,
const StdRegions::ConstFactorMap factors = StdRegions::NullConstFactorMap,
const StdRegions::VarCoeffMap varcoeffs = StdRegions::NullVarCoeffMap 
)

Definition at line 158 of file Expansion.cpp.

162{
163 MatrixKey mkey(mtype, DetShapeType(), *this, factors, varcoeffs);
164 return GetLocMatrix(mkey);
165}
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.

References Nektar::StdRegions::StdExpansion::DetShapeType(), and GetLocMatrix().

◆ GetMetricInfo()

const SpatialDomains::GeomFactorsSharedPtr & Nektar::LocalRegions::Expansion::GetMetricInfo ( ) const

◆ GetMF()

Array< OneD, NekDouble > Nektar::LocalRegions::Expansion::GetMF ( const int  dir,
const int  shapedim,
const StdRegions::VarCoeffMap varcoeffs 
)
protected

Definition at line 718 of file Expansion.cpp.

720{
721 int coordim = GetCoordim();
722
723 int nquad0, nquad1, nquad2;
724 int nqtot = 1;
725 switch (shapedim)
726 {
727 case 2:
728 {
729 nquad0 = m_base[0]->GetNumPoints();
730 nquad1 = m_base[1]->GetNumPoints();
731 nqtot = nquad0 * nquad1;
732 break;
733 }
734 case 3:
735 {
736 nquad0 = m_base[0]->GetNumPoints();
737 nquad1 = m_base[1]->GetNumPoints();
738 nquad2 = m_base[2]->GetNumPoints();
739 nqtot = nquad0 * nquad1 * nquad2;
740 break;
741 }
742 default:
743 break;
744 }
745
746 StdRegions::VarCoeffType MMFCoeffs[15] = {
755
756 Array<OneD, NekDouble> MF(coordim * nqtot);
757 Array<OneD, NekDouble> tmp(nqtot);
758 for (int k = 0; k < coordim; k++)
759 {
760 StdRegions::VarCoeffMap::const_iterator MFdir =
761 varcoeffs.find(MMFCoeffs[5 * dir + k]);
762 tmp = MFdir->second.GetValue();
763
764 Vmath::Vcopy(nqtot, &tmp[0], 1, &MF[k * nqtot], 1);
765 }
766
767 return MF;
768}
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition Vmath.hpp:825

References Nektar::StdRegions::eVarCoeffMF1Div, Nektar::StdRegions::eVarCoeffMF1Mag, Nektar::StdRegions::eVarCoeffMF1x, Nektar::StdRegions::eVarCoeffMF1y, Nektar::StdRegions::eVarCoeffMF1z, Nektar::StdRegions::eVarCoeffMF2Div, Nektar::StdRegions::eVarCoeffMF2Mag, Nektar::StdRegions::eVarCoeffMF2x, Nektar::StdRegions::eVarCoeffMF2y, Nektar::StdRegions::eVarCoeffMF2z, Nektar::StdRegions::eVarCoeffMF3Div, Nektar::StdRegions::eVarCoeffMF3Mag, Nektar::StdRegions::eVarCoeffMF3x, Nektar::StdRegions::eVarCoeffMF3y, Nektar::StdRegions::eVarCoeffMF3z, Nektar::StdRegions::StdExpansion::GetCoordim(), Nektar::StdRegions::StdExpansion::m_base, and Vmath::Vcopy().

Referenced by Nektar::LocalRegions::Expansion2D::AddHDGHelmholtzEdgeTerms(), Nektar::LocalRegions::Expansion3D::AddHDGHelmholtzFaceTerms(), Nektar::LocalRegions::Expansion2D::v_GenMatrix(), and Nektar::LocalRegions::Expansion3D::v_GenMatrix().

◆ GetMFDiv()

Array< OneD, NekDouble > Nektar::LocalRegions::Expansion::GetMFDiv ( const int  dir,
const StdRegions::VarCoeffMap varcoeffs 
)
protected

Definition at line 771 of file Expansion.cpp.

773{
774 int indxDiv = 3;
775
776 StdRegions::VarCoeffType MMFCoeffs[15] = {
785
786 StdRegions::VarCoeffMap::const_iterator MFdir =
787 varcoeffs.find(MMFCoeffs[5 * dir + indxDiv]);
788 Array<OneD, NekDouble> MFDiv = MFdir->second.GetValue();
789
790 return MFDiv;
791}

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, and Nektar::StdRegions::eVarCoeffMF3z.

Referenced by Nektar::LocalRegions::Expansion2D::AddHDGHelmholtzEdgeTerms(), Nektar::LocalRegions::Expansion3D::AddHDGHelmholtzFaceTerms(), Nektar::LocalRegions::Expansion2D::v_GenMatrix(), and Nektar::LocalRegions::Expansion3D::v_GenMatrix().

◆ GetMFMag()

Array< OneD, NekDouble > Nektar::LocalRegions::Expansion::GetMFMag ( const int  dir,
const StdRegions::VarCoeffMap varcoeffs 
)
protected

Definition at line 794 of file Expansion.cpp.

796{
797 int indxMag = 4;
798
799 StdRegions::VarCoeffType MMFCoeffs[15] = {
808
809 StdRegions::VarCoeffMap::const_iterator MFdir =
810 varcoeffs.find(MMFCoeffs[5 * dir + indxMag]);
811 Array<OneD, NekDouble> MFmag = MFdir->second.GetValue();
812
813 return MFmag;
814}

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, and Nektar::StdRegions::eVarCoeffMF3z.

Referenced by Nektar::LocalRegions::Expansion2D::AddHDGHelmholtzEdgeTerms(), Nektar::LocalRegions::Expansion3D::AddHDGHelmholtzFaceTerms(), Nektar::LocalRegions::Expansion2D::v_GenMatrix(), and Nektar::LocalRegions::Expansion3D::v_GenMatrix().

◆ GetPhysNormals()

const Array< OneD, const NekDouble > & Nektar::LocalRegions::Expansion::GetPhysNormals ( void  )
inline

Definition at line 228 of file Expansion.h.

229 {
230 return v_GetPhysNormals();
231 }
virtual const Array< OneD, const NekDouble > & v_GetPhysNormals()

References v_GetPhysNormals().

◆ GetRightAdjacentElementExp()

ExpansionSharedPtr Nektar::LocalRegions::Expansion::GetRightAdjacentElementExp ( ) const
inline

Definition at line 453 of file Expansion.h.

454{
455 ASSERTL1(m_elementLeft.lock().get(), "Right adjacent element not set.");
456
457 return m_elementRight.lock();
458}
ExpansionWeakPtr m_elementRight
Definition Expansion.h:284

References ASSERTL1, m_elementLeft, and m_elementRight.

◆ GetRightAdjacentElementTrace()

int Nektar::LocalRegions::Expansion::GetRightAdjacentElementTrace ( ) const
inline

Definition at line 465 of file Expansion.h.

466{
467 return m_elementTraceRight;
468}

References m_elementTraceRight.

◆ GetTraceExp()

ExpansionSharedPtr Nektar::LocalRegions::Expansion::GetTraceExp ( const int  traceid)
inline

Definition at line 420 of file Expansion.h.

421{
422 ASSERTL1(traceid < GetNtraces(), "Trace is out of range.");
423
424 ExpansionSharedPtr returnval;
425
426 if (m_traceExp.count(traceid))
427 {
428 // Use stored value
429 returnval = m_traceExp[traceid].lock();
430 }
431 else
432 {
433 // Generate trace exp
434 v_GenTraceExp(traceid, returnval);
435 }
436
437 return returnval;
438}
virtual void v_GenTraceExp(const int traceid, ExpansionSharedPtr &exp)
int GetNtraces() const
Returns the number of trace elements connected to this element.
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition Expansion.h:66

References ASSERTL1, Nektar::StdRegions::StdExpansion::GetNtraces(), m_traceExp, and v_GenTraceExp().

Referenced by Nektar::LocalRegions::Expansion3D::AddHDGHelmholtzFaceTerms(), Nektar::LocalRegions::Expansion2D::v_GenMatrix(), and Nektar::LocalRegions::Expansion3D::v_GenMatrix().

◆ GetTraceNormal()

const NormalVector & Nektar::LocalRegions::Expansion::GetTraceNormal ( const int  id)

◆ GetTraceOrient()

StdRegions::Orientation Nektar::LocalRegions::Expansion::GetTraceOrient ( int  trace)
inline

◆ GetTracePhysMap()

void Nektar::LocalRegions::Expansion::GetTracePhysMap ( const int  edge,
Array< OneD, int > &  outarray 
)
inline

Definition at line 209 of file Expansion.h.

210 {
211 v_GetTracePhysMap(edge, outarray);
212 }
virtual void v_GetTracePhysMap(const int edge, Array< OneD, int > &outarray)

References v_GetTracePhysMap().

Referenced by PhysBaseOnTraceMat(), PhysDerivBaseOnTraceMat(), StdDerivBaseOnTraceMat(), and Nektar::LocalRegions::Expansion3D::v_GetTracePhysVals().

◆ GetTracePhysVals()

void Nektar::LocalRegions::Expansion::GetTracePhysVals ( const int  trace,
const StdRegions::StdExpansionSharedPtr TraceExp,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
StdRegions::Orientation  orient = StdRegions::eNoOrientation 
)
inline

Definition at line 200 of file Expansion.h.

205 {
206 v_GetTracePhysVals(trace, TraceExp, inarray, outarray, orient);
207 }
virtual void v_GetTracePhysVals(const int trace, const StdRegions::StdExpansionSharedPtr &TraceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient)

References v_GetTracePhysVals().

Referenced by Nektar::LocalRegions::Expansion2D::v_GenMatrix(), and Nektar::LocalRegions::Expansion3D::v_GenMatrix().

◆ GetTraceQFactors()

void Nektar::LocalRegions::Expansion::GetTraceQFactors ( const int  trace,
Array< OneD, NekDouble > &  outarray 
)
inline

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

Definition at line 194 of file Expansion.h.

196 {
197 v_GetTraceQFactors(trace, outarray);
198 }
virtual void v_GetTraceQFactors(const int trace, Array< OneD, NekDouble > &outarray)

References v_GetTraceQFactors().

◆ NormalTraceDerivFactors()

void Nektar::LocalRegions::Expansion::NormalTraceDerivFactors ( Array< OneD, Array< OneD, NekDouble > > &  factors,
Array< OneD, Array< OneD, NekDouble > > &  d0factors,
Array< OneD, Array< OneD, NekDouble > > &  d1factors 
)

Definition at line 150 of file Expansion.cpp.

154{
155 return v_NormalTraceDerivFactors(factors, d0factors, d1factors);
156}
virtual void v_NormalTraceDerivFactors(Array< OneD, Array< OneD, NekDouble > > &factors, Array< OneD, Array< OneD, NekDouble > > &d0factors, Array< OneD, Array< OneD, NekDouble > > &d1factors)

References v_NormalTraceDerivFactors().

◆ PhysBaseOnTraceMat()

void Nektar::LocalRegions::Expansion::PhysBaseOnTraceMat ( const int  traceid,
DNekMatSharedPtr BdataMat 
)

Definition at line 586 of file Expansion.cpp.

588{
589 int nquad = GetTotPoints();
590
591 Array<OneD, NekDouble> coeffs(m_ncoeffs);
592 Array<OneD, NekDouble> phys(nquad);
593
594 Array<OneD, int> tracePhysIds;
595 GetTracePhysMap(traceid, tracePhysIds);
596
597 int nTracePts = GetTraceNumPoints(traceid);
598
599 BdataMat =
601
602 for (int i = 0; i < m_ncoeffs; ++i)
603 {
604 Vmath::Zero(m_ncoeffs, coeffs, 1);
605 coeffs[i] = 1.0;
606 BwdTrans(coeffs, phys);
607
608 for (int k = 0; k < nTracePts; ++k)
609 {
610 (*BdataMat)(i, k) = phys[tracePhysIds[k]];
611 }
612 }
613}
void GetTracePhysMap(const int edge, Array< OneD, int > &outarray)
Definition Expansion.h:209
int GetTraceNumPoints(const int i) const
This function returns the number of quadrature points belonging to the i-th trace.
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Backward transformation from coefficient space to physical space.
void Zero(int n, T *x, const int incx)
Zero vector.
Definition Vmath.hpp:273

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::StdRegions::StdExpansion::BwdTrans(), Nektar::StdRegions::StdExpansion::GetTotPoints(), Nektar::StdRegions::StdExpansion::GetTraceNumPoints(), GetTracePhysMap(), Nektar::StdRegions::StdExpansion::m_ncoeffs, and Vmath::Zero().

◆ PhysDerivBaseOnTraceMat()

void Nektar::LocalRegions::Expansion::PhysDerivBaseOnTraceMat ( const int  traceid,
Array< OneD, DNekMatSharedPtr > &  DerivMat 
)

Definition at line 542 of file Expansion.cpp.

544{
545 int nquad = GetTotPoints();
546 int ndir = m_base.size();
547
548 Array<OneD, NekDouble> coeffs(m_ncoeffs);
549 Array<OneD, NekDouble> phys(nquad);
550
551 Array<OneD, int> tracePhysIds;
552 GetTracePhysMap(traceid, tracePhysIds);
553
554 int nTracePts = GetTraceNumPoints(traceid);
555
556 // initialise array to null so can call for
557 // differnt dimensions
558 Array<OneD, Array<OneD, NekDouble>> Deriv(3, NullNekDouble1DArray);
559 DerivMat = Array<OneD, DNekMatSharedPtr>(ndir);
560 for (int d = 0; d < ndir; ++d)
561 {
562 Deriv[d] = Array<OneD, NekDouble>(nquad);
564 nTracePts, 0.0);
565 }
566
567 for (int i = 0; i < m_ncoeffs; ++i)
568 {
569 Vmath::Zero(m_ncoeffs, coeffs, 1);
570 coeffs[i] = 1.0;
571 BwdTrans(coeffs, phys);
572
573 // dphi_i/d\xi_1, dphi_i/d\xi_2 dphi_i/d\xi_3
574 PhysDeriv(phys, Deriv[0], Deriv[1], Deriv[2]);
575
576 for (int k = 0; k < nTracePts; ++k)
577 {
578 for (int d = 0; d < ndir; ++d)
579 {
580 (*DerivMat[d])(i, k) = Deriv[d][tracePhysIds[k]];
581 }
582 }
583 }
584}
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)
std::vector< double > d(NPUPPER *NPUPPER)
static Array< OneD, NekDouble > NullNekDouble1DArray

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::StdRegions::StdExpansion::BwdTrans(), Nektar::StdRegions::StdExpansion::GetTotPoints(), Nektar::StdRegions::StdExpansion::GetTraceNumPoints(), GetTracePhysMap(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::NullNekDouble1DArray, Nektar::StdRegions::StdExpansion::PhysDeriv(), and Vmath::Zero().

◆ ReOrientTracePhysMap()

void Nektar::LocalRegions::Expansion::ReOrientTracePhysMap ( const StdRegions::Orientation  orient,
Array< OneD, int > &  idmap,
const int  nq0,
const int  nq1 
)
inline

Definition at line 214 of file Expansion.h.

217 {
218 v_ReOrientTracePhysMap(orient, idmap, nq0, nq1);
219 }
virtual void v_ReOrientTracePhysMap(const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1=-1)

References v_ReOrientTracePhysMap().

◆ Reset()

void Nektar::LocalRegions::Expansion::Reset ( )

Definition at line 172 of file Expansion.cpp.

173{
174 // Clear metrics
175 m_metrics.clear();
176
177 // Regenerate geometry factors
179}

References Nektar::SpatialDomains::Geometry::GetGeomFactors(), m_geom, m_metricinfo, and m_metrics.

Referenced by export_Expansion().

◆ SetAdjacentElementExp()

void Nektar::LocalRegions::Expansion::SetAdjacentElementExp ( int  traceid,
ExpansionSharedPtr e 
)
inline

Definition at line 470 of file Expansion.h.

472{
473 if (m_elementLeft.lock().get())
474 {
475 m_elementRight = exp;
476 m_elementTraceRight = traceid;
477 }
478 else
479 {
480 m_elementLeft = exp;
481 m_elementTraceLeft = traceid;
482 }
483}

References m_elementLeft, m_elementRight, m_elementTraceLeft, and m_elementTraceRight.

◆ SetCoeffsToOrientation()

void Nektar::LocalRegions::Expansion::SetCoeffsToOrientation ( StdRegions::Orientation  dir,
Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
inline

Definition at line 173 of file Expansion.h.

176 {
177 v_SetCoeffsToOrientation(dir, inarray, outarray);
178 }
void v_SetCoeffsToOrientation(StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override

References v_SetCoeffsToOrientation().

◆ SetPhysNormals()

void Nektar::LocalRegions::Expansion::SetPhysNormals ( Array< OneD, const NekDouble > &  normal)
inline

Definition at line 233 of file Expansion.h.

234 {
235 v_SetPhysNormals(normal);
236 }
virtual void v_SetPhysNormals(Array< OneD, const NekDouble > &normal)

References v_SetPhysNormals().

◆ SetTraceExp()

void Nektar::LocalRegions::Expansion::SetTraceExp ( const int  traceid,
ExpansionSharedPtr f 
)
inline

Definition at line 440 of file Expansion.h.

441{
442 ASSERTL1(traceid < GetNtraces(), "Trace out of range.");
443
444 m_traceExp[traceid] = exp;
445}

References ASSERTL1, Nektar::StdRegions::StdExpansion::GetNtraces(), and m_traceExp.

◆ SetUpPhysNormals()

void Nektar::LocalRegions::Expansion::SetUpPhysNormals ( const int  trace)
inline

Definition at line 238 of file Expansion.h.

239 {
240 v_SetUpPhysNormals(trace);
241 }
virtual void v_SetUpPhysNormals(const int id)

References v_SetUpPhysNormals().

◆ StdDerivBaseOnTraceMat()

void Nektar::LocalRegions::Expansion::StdDerivBaseOnTraceMat ( Array< OneD, DNekMatSharedPtr > &  DerivMat)

Definition at line 486 of file Expansion.cpp.

487{
488 int nquad = GetTotPoints();
489 int ntraces = GetNtraces();
490 int ndir = m_base.size();
491
492 Array<OneD, NekDouble> coeffs(m_ncoeffs);
493 Array<OneD, NekDouble> phys(nquad);
494
495 Array<OneD, Array<OneD, int>> traceids(ntraces);
496
497 int tottracepts = 0;
498 for (int i = 0; i < ntraces; ++i)
499 {
500 GetTracePhysMap(i, traceids[i]);
501 tottracepts += GetTraceNumPoints(i);
502 }
503
504 // initialise array to null so can call for
505 // differnt dimensions
506 Array<OneD, Array<OneD, NekDouble>> Deriv(3, NullNekDouble1DArray);
507
508 DerivMat = Array<OneD, DNekMatSharedPtr>(ndir);
509
510 for (int i = 0; i < ndir; ++i)
511 {
512 Deriv[i] = Array<OneD, NekDouble>(nquad);
513 DerivMat[i] =
515 }
516
517 for (int i = 0; i < m_ncoeffs; ++i)
518 {
519 Vmath::Zero(m_ncoeffs, coeffs, 1);
520 coeffs[i] = 1.0;
521 BwdTrans(coeffs, phys);
522
523 // dphi_i/d\xi_1, dphi_i/d\xi_2 dphi_i/d\xi_3
524 StdPhysDeriv(phys, Deriv[0], Deriv[1], Deriv[2]);
525
526 int cnt = 0;
527 for (int j = 0; j < ntraces; ++j)
528 {
529 int nTracePts = GetTraceNumPoints(j);
530 for (int k = 0; k < nTracePts; ++k)
531 {
532 for (int d = 0; d < ndir; ++d)
533 {
534 (*DerivMat[d])(i, cnt + k) = Deriv[d][traceids[j][k]];
535 }
536 }
537 cnt += nTracePts;
538 }
539 }
540}
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)

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::StdRegions::StdExpansion::BwdTrans(), Nektar::StdRegions::StdExpansion::GetNtraces(), Nektar::StdRegions::StdExpansion::GetTotPoints(), Nektar::StdRegions::StdExpansion::GetTraceNumPoints(), GetTracePhysMap(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::NullNekDouble1DArray, Nektar::StdRegions::StdExpansion::StdPhysDeriv(), and Vmath::Zero().

◆ TraceNormLen()

void Nektar::LocalRegions::Expansion::TraceNormLen ( const int  traceid,
NekDouble h,
NekDouble p 
)
inline

Definition at line 250 of file Expansion.h.

251 {
252 v_TraceNormLen(traceid, h, p);
253 }
virtual void v_TraceNormLen(const int traceid, NekDouble &h, NekDouble &p)

References v_TraceNormLen().

Referenced by Nektar::LocalRegions::Expansion2D::v_GenMatrix(), and Nektar::LocalRegions::Expansion3D::v_GenMatrix().

◆ v_AddEdgeNormBoundaryInt() [1/2]

void Nektar::LocalRegions::Expansion::v_AddEdgeNormBoundaryInt ( const int  edge,
const std::shared_ptr< Expansion > &  EdgeExp,
const Array< OneD, const NekDouble > &  Fn,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Definition at line 851 of file Expansion.cpp.

856{
857 NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
858}

References Nektar::ErrorUtil::efatal, and NEKERROR.

◆ v_AddEdgeNormBoundaryInt() [2/2]

void Nektar::LocalRegions::Expansion::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 
)
protectedvirtual

Definition at line 841 of file Expansion.cpp.

847{
848 NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
849}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by AddEdgeNormBoundaryInt(), and AddEdgeNormBoundaryInt().

◆ v_AddFaceNormBoundaryInt()

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

Definition at line 860 of file Expansion.cpp.

865{
866 NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
867}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by AddFaceNormBoundaryInt().

◆ v_AddRobinMassMatrix()

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

Reimplemented in Nektar::LocalRegions::Expansion2D, Nektar::LocalRegions::Expansion3D, and Nektar::LocalRegions::Expansion1D.

Definition at line 975 of file Expansion.cpp.

979{
980 NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
981}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by AddRobinMassMatrix().

◆ v_AddRobinTraceContribution()

void Nektar::LocalRegions::Expansion::v_AddRobinTraceContribution ( const int  traceid,
const Array< OneD, const NekDouble > &  primCoeffs,
const Array< OneD, NekDouble > &  incoeffs,
Array< OneD, NekDouble > &  coeffs 
)
protectedvirtual

Reimplemented in Nektar::LocalRegions::Expansion2D, and Nektar::LocalRegions::Expansion1D.

Definition at line 983 of file Expansion.cpp.

988{
989 NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
990}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by AddRobinTraceContribution().

◆ v_AlignVectorToCollapsedDir()

void Nektar::LocalRegions::Expansion::v_AlignVectorToCollapsedDir ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protectedvirtual

◆ v_BuildTransformationMatrix()

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

Reimplemented in Nektar::LocalRegions::Expansion3D.

Definition at line 816 of file Expansion.cpp.

819{
820 NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
822}
static DNekMatSharedPtr NullDNekMatSharedPtr

References Nektar::ErrorUtil::efatal, NEKERROR, and Nektar::NullDNekMatSharedPtr.

Referenced by BuildTransformationMatrix().

◆ v_BuildVertexMatrix()

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

Reimplemented in Nektar::LocalRegions::Expansion2D, and Nektar::LocalRegions::Expansion3D.

Definition at line 824 of file Expansion.cpp.

826{
827 NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
829}

References Nektar::ErrorUtil::efatal, NEKERROR, and Nektar::NullDNekMatSharedPtr.

Referenced by BuildVertexMatrix().

◆ v_ComputeLaplacianMetric()

virtual void Nektar::LocalRegions::Expansion::v_ComputeLaplacianMetric ( )
inlineprotectedvirtual

◆ v_ComputeTraceNormal()

void Nektar::LocalRegions::Expansion::v_ComputeTraceNormal ( const int  id)
protectedvirtual

◆ v_DGDeriv()

void Nektar::LocalRegions::Expansion::v_DGDeriv ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, ExpansionSharedPtr > &  EdgeExp,
Array< OneD, Array< OneD, NekDouble > > &  coeffs,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented in Nektar::LocalRegions::Expansion2D, and Nektar::LocalRegions::Expansion3D.

Definition at line 869 of file Expansion.cpp.

875{
876 NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
877}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by DGDeriv().

◆ v_DivideByQuadratureMetric()

void Nektar::LocalRegions::Expansion::v_DivideByQuadratureMetric ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Definition at line 438 of file Expansion.cpp.

441{
442 const int nqtot = GetTotPoints();
443
444 if (m_metrics.count(eMetricQuadrature) == 0)
445 {
447 }
448
449 Vmath::Vdiv(nqtot, inarray, 1, m_metrics[eMetricQuadrature], 1, outarray,
450 1);
451 // remove NaN or Inf values and set to zero
452 for (int i = 0; i < nqtot; ++i)
453 {
454 if (std::isnan(outarray[i]) || std::isinf(outarray[i]))
455 {
456 outarray[i] = 0.0;
457 }
458 }
459}
void Vdiv(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.hpp:126

References ComputeQuadratureMetric(), Nektar::LocalRegions::eMetricQuadrature, Nektar::StdRegions::StdExpansion::GetTotPoints(), m_metrics, and Vmath::Vdiv().

Referenced by DivideByQuadratureMetric().

◆ v_DropLocMatrix()

void Nektar::LocalRegions::Expansion::v_DropLocMatrix ( const LocalRegions::MatrixKey mkey)
protectedvirtual

◆ v_ExtractDataToCoeffs()

void Nektar::LocalRegions::Expansion::v_ExtractDataToCoeffs ( const NekDouble data,
const std::vector< unsigned int > &  nummodes,
const int  nmodes_offset,
NekDouble coeffs,
std::vector< LibUtilities::BasisType > &  fromType 
)
protectedvirtual

◆ v_GenTraceExp()

void Nektar::LocalRegions::Expansion::v_GenTraceExp ( const int  traceid,
ExpansionSharedPtr exp 
)
protectedvirtual

Reimplemented in Nektar::LocalRegions::Expansion2D, and Nektar::LocalRegions::Expansion3D.

Definition at line 946 of file Expansion.cpp.

948{
950 "Method does not exist for this shape or library");
951}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetTraceExp().

◆ v_GetCoordim()

int Nektar::LocalRegions::Expansion::v_GetCoordim ( ) const
inlineoverrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 321 of file Expansion.h.

322 {
323 return m_geom->GetCoordim();
324 }
int GetCoordim() const
Return the coordinate dimension of this object (i.e. the dimension of the space in which this object ...
Definition Geometry.h:279

References Nektar::SpatialDomains::Geometry::GetCoordim(), and m_geom.

◆ v_GetCoords()

void Nektar::LocalRegions::Expansion::v_GetCoords ( Array< OneD, NekDouble > &  coords_1,
Array< OneD, NekDouble > &  coords_2,
Array< OneD, NekDouble > &  coords_3 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::PointExp, Nektar::LocalRegions::NodalTriExp, Nektar::LocalRegions::HexExp, Nektar::LocalRegions::PrismExp, Nektar::LocalRegions::PyrExp, Nektar::LocalRegions::QuadExp, Nektar::LocalRegions::SegExp, Nektar::LocalRegions::TetExp, and Nektar::LocalRegions::TriExp.

Definition at line 615 of file Expansion.cpp.

618{
619 ASSERTL1(m_geom, "m_geom not defined");
620
621 // get physical points defined in Geom
622 m_geom->FillGeom();
623
624 const int expDim = m_base.size();
625 int nqGeom = 1;
626 bool doCopy = true;
627
628 Array<OneD, LibUtilities::BasisSharedPtr> CBasis(expDim);
629 Array<OneD, Array<OneD, NekDouble>> tmp(3);
630
631 for (int i = 0; i < expDim; ++i)
632 {
633 CBasis[i] = m_geom->GetXmap()->GetBasis(i);
634 nqGeom *= CBasis[i]->GetNumPoints();
635 doCopy = doCopy &&
636 m_base[i]->GetBasisKey().SamePoints(CBasis[i]->GetBasisKey());
637 }
638
639 tmp[0] = coords_0;
640 tmp[1] = coords_1;
641 tmp[2] = coords_2;
642
643 if (doCopy)
644 {
645 for (int i = 0; i < m_geom->GetCoordim(); ++i)
646 {
647 m_geom->GetXmap()->BwdTrans(m_geom->GetCoeffs(i), tmp[i]);
648 }
649 }
650 else
651 {
652 for (int i = 0; i < m_geom->GetCoordim(); ++i)
653 {
654 Array<OneD, NekDouble> tmpGeom(nqGeom);
655 m_geom->GetXmap()->BwdTrans(m_geom->GetCoeffs(i), tmpGeom);
656
657 switch (expDim)
658 {
659 case 1:
660 {
662 CBasis[0]->GetPointsKey(), &tmpGeom[0],
663 m_base[0]->GetPointsKey(), &tmp[i][0]);
664 break;
665 }
666 case 2:
667 {
669 CBasis[0]->GetPointsKey(), CBasis[1]->GetPointsKey(),
670 &tmpGeom[0], m_base[0]->GetPointsKey(),
671 m_base[1]->GetPointsKey(), &tmp[i][0]);
672 break;
673 }
674 case 3:
675 {
677 CBasis[0]->GetPointsKey(), CBasis[1]->GetPointsKey(),
678 CBasis[2]->GetPointsKey(), &tmpGeom[0],
679 m_base[0]->GetPointsKey(), m_base[1]->GetPointsKey(),
680 m_base[2]->GetPointsKey(), &tmp[i][0]);
681 break;
682 }
683 }
684 }
685 }
686}
const Array< OneD, const NekDouble > & GetCoeffs(const int i) const
Return the coefficients of the transformation Geometry::m_xmap in coordinate direction i.
Definition Geometry.h:448
void FillGeom()
Populate the coordinate mapping Geometry::m_coeffs information from any children geometry elements.
Definition Geometry.h:460
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
Definition Geometry.h:439
void Interp1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
this function interpolates a 1D function evaluated at the quadrature points of the basis fbasis0 to ...
Definition Interp.cpp:47
void Interp3D(const BasisKey &fbasis0, const BasisKey &fbasis1, const BasisKey &fbasis2, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, const BasisKey &tbasis2, Array< OneD, NekDouble > &to)
this function interpolates a 3D function evaluated at the quadrature points of the 3D basis,...
Definition Interp.cpp:162
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:101

References ASSERTL1, Nektar::SpatialDomains::Geometry::FillGeom(), Nektar::SpatialDomains::Geometry::GetCoeffs(), Nektar::SpatialDomains::Geometry::GetCoordim(), Nektar::SpatialDomains::Geometry::GetXmap(), Nektar::LibUtilities::Interp1D(), Nektar::LibUtilities::Interp2D(), Nektar::LibUtilities::Interp3D(), Nektar::StdRegions::StdExpansion::m_base, and m_geom.

Referenced by Nektar::LocalRegions::NodalTriExp::v_GetCoords(), Nektar::LocalRegions::HexExp::v_GetCoords(), Nektar::LocalRegions::PrismExp::v_GetCoords(), Nektar::LocalRegions::PyrExp::v_GetCoords(), Nektar::LocalRegions::QuadExp::v_GetCoords(), Nektar::LocalRegions::SegExp::v_GetCoords(), Nektar::LocalRegions::TetExp::v_GetCoords(), and Nektar::LocalRegions::TriExp::v_GetCoords().

◆ v_GetLocMatrix()

DNekScalMatSharedPtr Nektar::LocalRegions::Expansion::v_GetLocMatrix ( const LocalRegions::MatrixKey mkey)
protectedvirtual

◆ v_GetPhysNormals()

const Array< OneD, const NekDouble > & Nektar::LocalRegions::Expansion::v_GetPhysNormals ( void  )
protectedvirtual

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 958 of file Expansion.cpp.

959{
960 NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
962}

References Nektar::ErrorUtil::efatal, NEKERROR, and Nektar::NullNekDouble1DArray.

Referenced by GetPhysNormals().

◆ v_GetTraceOrient()

StdRegions::Orientation Nektar::LocalRegions::Expansion::v_GetTraceOrient ( int  trace)
protectedvirtual

◆ v_GetTracePhysMap()

void Nektar::LocalRegions::Expansion::v_GetTracePhysMap ( const int  edge,
Array< OneD, int > &  outarray 
)
protectedvirtual

◆ v_GetTracePhysVals()

void Nektar::LocalRegions::Expansion::v_GetTracePhysVals ( const int  trace,
const StdRegions::StdExpansionSharedPtr TraceExp,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
StdRegions::Orientation  orient 
)
protectedvirtual

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

Definition at line 919 of file Expansion.cpp.

925{
927 "Method does not exist for this shape or library");
928}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetTracePhysVals().

◆ v_GetTraceQFactors()

void Nektar::LocalRegions::Expansion::v_GetTraceQFactors ( const int  trace,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

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

Definition at line 911 of file Expansion.cpp.

914{
916 "Method does not exist for this shape or library");
917}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetTraceQFactors().

◆ v_MultiplyByQuadratureMetric()

void Nektar::LocalRegions::Expansion::v_MultiplyByQuadratureMetric ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 423 of file Expansion.cpp.

426{
427 const int nqtot = GetTotPoints();
428
429 if (m_metrics.count(eMetricQuadrature) == 0)
430 {
432 }
433
434 Vmath::Vmul(nqtot, m_metrics[eMetricQuadrature], 1, inarray, 1, outarray,
435 1);
436}
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.hpp:72

References ComputeQuadratureMetric(), Nektar::LocalRegions::eMetricQuadrature, Nektar::StdRegions::StdExpansion::GetTotPoints(), m_metrics, and Vmath::Vmul().

◆ v_NormalTraceDerivFactors()

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

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

Definition at line 888 of file Expansion.cpp.

892{
894 "This function is only valid for "
895 "shape expansion in LocalRegions, not parant class");
896}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by NormalTraceDerivFactors().

◆ v_ReOrientTracePhysMap()

void Nektar::LocalRegions::Expansion::v_ReOrientTracePhysMap ( const StdRegions::Orientation  orient,
Array< OneD, int > &  idmap,
const int  nq0,
const int  nq1 = -1 
)
protectedvirtual

Reimplemented in Nektar::LocalRegions::Expansion1D, Nektar::LocalRegions::Expansion2D, and Nektar::LocalRegions::Expansion3D.

Definition at line 937 of file Expansion.cpp.

941{
943 "Method does not exist for this shape or library");
944}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by ReOrientTracePhysMap().

◆ v_SetCoeffsToOrientation()

void Nektar::LocalRegions::Expansion::v_SetCoeffsToOrientation ( StdRegions::Orientation  dir,
Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 903 of file Expansion.cpp.

907{
908 NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
909}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by SetCoeffsToOrientation().

◆ v_SetPhysNormals()

void Nektar::LocalRegions::Expansion::v_SetPhysNormals ( Array< OneD, const NekDouble > &  normal)
protectedvirtual

Definition at line 964 of file Expansion.cpp.

966{
967 NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
968}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by SetPhysNormals().

◆ v_SetUpPhysNormals()

void Nektar::LocalRegions::Expansion::v_SetUpPhysNormals ( const int  id)
protectedvirtual

Reimplemented in Nektar::LocalRegions::Expansion2D.

Definition at line 970 of file Expansion.cpp.

971{
972 NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
973}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by SetUpPhysNormals().

◆ v_TraceNormLen()

void Nektar::LocalRegions::Expansion::v_TraceNormLen ( const int  traceid,
NekDouble h,
NekDouble p 
)
protectedvirtual

Reimplemented in Nektar::LocalRegions::Expansion1D, Nektar::LocalRegions::Expansion2D, and Nektar::LocalRegions::Expansion3D.

Definition at line 992 of file Expansion.cpp.

995{
996 NEKERROR(ErrorUtil::efatal, "This method has not been defined");
997}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by TraceNormLen().

◆ v_VectorFlux()

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

Reimplemented in Nektar::LocalRegions::Expansion1D, and Nektar::LocalRegions::Expansion2D.

Definition at line 879 of file Expansion.cpp.

881{
883 "This function is only valid for "
884 "shape expansion in LocalRegions, not parant class");
885 return 0.0;
886}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by VectorFlux().

◆ VectorFlux()

NekDouble Nektar::LocalRegions::Expansion::VectorFlux ( const Array< OneD, Array< OneD, NekDouble > > &  vec)

Definition at line 145 of file Expansion.cpp.

146{
147 return v_VectorFlux(vec);
148}
virtual NekDouble v_VectorFlux(const Array< OneD, Array< OneD, NekDouble > > &vec)

References v_VectorFlux().

Member Data Documentation

◆ m_elementLeft

ExpansionWeakPtr Nektar::LocalRegions::Expansion::m_elementLeft
protected

◆ m_elementRight

ExpansionWeakPtr Nektar::LocalRegions::Expansion::m_elementRight
protected

Definition at line 284 of file Expansion.h.

Referenced by GetRightAdjacentElementExp(), and SetAdjacentElementExp().

◆ m_elementTraceLeft

int Nektar::LocalRegions::Expansion::m_elementTraceLeft = -1
protected

Definition at line 285 of file Expansion.h.

Referenced by GetLeftAdjacentElementTrace(), and SetAdjacentElementExp().

◆ m_elementTraceRight

int Nektar::LocalRegions::Expansion::m_elementTraceRight = -1
protected

Definition at line 286 of file Expansion.h.

Referenced by GetRightAdjacentElementTrace(), and SetAdjacentElementExp().

◆ m_elmtBndNormDirElmtLen

std::map<int, Array<OneD, NekDouble> > Nektar::LocalRegions::Expansion::m_elmtBndNormDirElmtLen
protected

◆ m_geom

SpatialDomains::Geometry* Nektar::LocalRegions::Expansion::m_geom
protected

Definition at line 279 of file Expansion.h.

Referenced by Nektar::LocalRegions::SegExp::CreateMatrix(), Expansion(), Nektar::LocalRegions::PointExp::GetCoords(), GetGeom(), Nektar::LocalRegions::PointExp::GetGeom(), Nektar::LocalRegions::Expansion0D::GetGeom0D(), Nektar::LocalRegions::Expansion2D::GetGeom2D(), Nektar::LocalRegions::Expansion3D::GetGeom3D(), Reset(), Nektar::LocalRegions::Expansion2D::v_AddEdgeNormBoundaryInt(), Nektar::LocalRegions::NodalTriExp::v_AlignVectorToCollapsedDir(), Nektar::LocalRegions::QuadExp::v_AlignVectorToCollapsedDir(), Nektar::LocalRegions::TriExp::v_AlignVectorToCollapsedDir(), Nektar::LocalRegions::Expansion2D::v_GenTraceExp(), Nektar::LocalRegions::Expansion3D::v_GenTraceExp(), Nektar::LocalRegions::NodalTriExp::v_GetCoord(), Nektar::LocalRegions::HexExp::v_GetCoord(), Nektar::LocalRegions::PrismExp::v_GetCoord(), Nektar::LocalRegions::PyrExp::v_GetCoord(), Nektar::LocalRegions::QuadExp::v_GetCoord(), Nektar::LocalRegions::SegExp::v_GetCoord(), Nektar::LocalRegions::TetExp::v_GetCoord(), Nektar::LocalRegions::TriExp::v_GetCoord(), v_GetCoordim(), Nektar::LocalRegions::PointExp::v_GetCoords(), v_GetCoords(), Nektar::LocalRegions::PrismExp::v_GetSimplexEquiSpacedConnectivity(), Nektar::LocalRegions::QuadExp::v_GetTraceOrient(), Nektar::LocalRegions::Expansion3D::v_GetTraceOrient(), Nektar::LocalRegions::SegExp::v_HelmholtzMatrixOp(), Nektar::LocalRegions::SegExp::v_IProductWRTDerivBase(), Nektar::LocalRegions::QuadExp::v_IProductWRTDerivBase_SumFac(), Nektar::LocalRegions::SegExp::v_LaplacianMatrixOp(), Nektar::LocalRegions::SegExp::v_PhysDeriv_n(), Nektar::LocalRegions::SegExp::v_PhysDeriv_s(), Nektar::LocalRegions::QuadExp::v_PhysDirectionalDeriv(), Nektar::LocalRegions::TriExp::v_PhysDirectionalDeriv(), Nektar::LocalRegions::HexExp::v_PhysEvalFirstDeriv(), Nektar::LocalRegions::PrismExp::v_PhysEvalFirstDeriv(), Nektar::LocalRegions::PyrExp::v_PhysEvalFirstDeriv(), Nektar::LocalRegions::QuadExp::v_PhysEvalFirstDeriv(), Nektar::LocalRegions::SegExp::v_PhysEvalFirstDeriv(), Nektar::LocalRegions::TetExp::v_PhysEvalFirstDeriv(), Nektar::LocalRegions::TriExp::v_PhysEvalFirstDeriv(), Nektar::LocalRegions::SegExp::v_PhysEvalFirstSecondDeriv(), Nektar::LocalRegions::NodalTriExp::v_PhysEvaluate(), Nektar::LocalRegions::PrismExp::v_PhysEvaluate(), Nektar::LocalRegions::PyrExp::v_PhysEvaluate(), Nektar::LocalRegions::QuadExp::v_PhysEvaluate(), Nektar::LocalRegions::SegExp::v_PhysEvaluate(), Nektar::LocalRegions::TriExp::v_PhysEvaluate(), Nektar::LocalRegions::HexExp::v_PhysEvaluate(), and Nektar::LocalRegions::TetExp::v_PhysEvaluate().

◆ m_indexMapManager

LibUtilities::NekManager<IndexMapKey, IndexMapValues, IndexMapKey::opLess> Nektar::LocalRegions::Expansion::m_indexMapManager
protected

Definition at line 276 of file Expansion.h.

Referenced by GetIndexMap().

◆ m_metricinfo

SpatialDomains::GeomFactorsSharedPtr Nektar::LocalRegions::Expansion::m_metricinfo
protected

Definition at line 280 of file Expansion.h.

Referenced by ComputeGmatcdotMF(), ComputeQuadratureMetric(), Nektar::LocalRegions::Expansion2D::CreateMatrix(), Nektar::LocalRegions::Expansion3D::CreateMatrix(), Nektar::LocalRegions::SegExp::CreateMatrix(), CreateStaticCondMatrix(), Expansion(), GetMetricInfo(), Reset(), Nektar::LocalRegions::HexExp::v_AlignVectorToCollapsedDir(), Nektar::LocalRegions::NodalTriExp::v_AlignVectorToCollapsedDir(), Nektar::LocalRegions::PrismExp::v_AlignVectorToCollapsedDir(), Nektar::LocalRegions::PyrExp::v_AlignVectorToCollapsedDir(), Nektar::LocalRegions::QuadExp::v_AlignVectorToCollapsedDir(), Nektar::LocalRegions::TetExp::v_AlignVectorToCollapsedDir(), Nektar::LocalRegions::TriExp::v_AlignVectorToCollapsedDir(), Nektar::LocalRegions::HexExp::v_ComputeLaplacianMetric(), Nektar::LocalRegions::PyrExp::v_ComputeLaplacianMetric(), Nektar::LocalRegions::QuadExp::v_ComputeLaplacianMetric(), Nektar::LocalRegions::TriExp::v_ComputeLaplacianMetric(), Nektar::LocalRegions::QuadExp::v_GetTraceQFactors(), Nektar::LocalRegions::SegExp::v_HelmholtzMatrixOp(), Nektar::LocalRegions::HexExp::v_Integral(), Nektar::LocalRegions::NodalTriExp::v_Integral(), Nektar::LocalRegions::PrismExp::v_Integral(), Nektar::LocalRegions::PyrExp::v_Integral(), Nektar::LocalRegions::QuadExp::v_Integral(), Nektar::LocalRegions::SegExp::v_Integral(), Nektar::LocalRegions::TetExp::v_Integral(), Nektar::LocalRegions::TriExp::v_Integral(), Nektar::LocalRegions::SegExp::v_IProductWRTBase(), Nektar::LocalRegions::SegExp::v_IProductWRTDerivBase(), Nektar::LocalRegions::HexExp::v_IProductWRTDirectionalDerivBase_SumFac(), Nektar::LocalRegions::TriExp::v_IProductWRTDirectionalDerivBase_SumFac(), Nektar::LocalRegions::SegExp::v_LaplacianMatrixOp(), Nektar::LocalRegions::PrismExp::v_LaplacianMatrixOp_MatFree_Kernel(), Nektar::LocalRegions::PrismExp::v_NormalTraceDerivFactors(), Nektar::LocalRegions::Expansion1D::v_NormalTraceDerivFactors(), Nektar::LocalRegions::HexExp::v_NormalTraceDerivFactors(), Nektar::LocalRegions::QuadExp::v_NormalTraceDerivFactors(), Nektar::LocalRegions::TriExp::v_NormalTraceDerivFactors(), Nektar::LocalRegions::QuadExp::v_NormVectorIProductWRTBase(), Nektar::LocalRegions::TriExp::v_NormVectorIProductWRTBase(), Nektar::LocalRegions::HexExp::v_PhysDeriv(), Nektar::LocalRegions::PrismExp::v_PhysDeriv(), Nektar::LocalRegions::PyrExp::v_PhysDeriv(), Nektar::LocalRegions::TetExp::v_PhysDeriv(), Nektar::LocalRegions::NodalTriExp::v_PhysDeriv(), Nektar::LocalRegions::QuadExp::v_PhysDeriv(), Nektar::LocalRegions::TriExp::v_PhysDeriv(), Nektar::LocalRegions::SegExp::v_PhysDeriv(), Nektar::LocalRegions::SegExp::v_PhysDeriv_n(), Nektar::LocalRegions::SegExp::v_PhysDeriv_s(), Nektar::LocalRegions::HexExp::v_PhysDirectionalDeriv(), Nektar::LocalRegions::QuadExp::v_PhysDirectionalDeriv(), Nektar::LocalRegions::TriExp::v_PhysDirectionalDeriv(), Nektar::LocalRegions::HexExp::v_SVVLaplacianFilter(), Nektar::LocalRegions::PrismExp::v_SVVLaplacianFilter(), Nektar::LocalRegions::PyrExp::v_SVVLaplacianFilter(), Nektar::LocalRegions::QuadExp::v_SVVLaplacianFilter(), Nektar::LocalRegions::TetExp::v_SVVLaplacianFilter(), and Nektar::LocalRegions::TriExp::v_SVVLaplacianFilter().

◆ m_metrics

MetricMap Nektar::LocalRegions::Expansion::m_metrics
protected

◆ m_traceExp

std::map<int, ExpansionWeakPtr> Nektar::LocalRegions::Expansion::m_traceExp
protected

◆ m_traceNormals

std::map<int, NormalVector> Nektar::LocalRegions::Expansion::m_traceNormals
protected