Nektar++
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::GeometrySharedPtr 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::GeometrySharedPtr GetGeom () const
 
void Reset ()
 
IndexMapValuesSharedPtr CreateIndexMap (const IndexMapKey &ikey)
 
DNekScalBlkMatSharedPtr CreateStaticCondMatrix (const MatrixKey &mkey)
 
const SpatialDomains::GeomFactorsSharedPtrGetMetricInfo () const
 
DNekMatSharedPtr BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
 
DNekMatSharedPtr BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd)
 
void ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType)
 
void AddEdgeNormBoundaryInt (const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
void AddEdgeNormBoundaryInt (const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
void AddFaceNormBoundaryInt (const int face, const std::shared_ptr< Expansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
void DGDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &coeffs, Array< OneD, NekDouble > &outarray)
 
NekDouble VectorFlux (const Array< OneD, Array< OneD, NekDouble > > &vec)
 
void NormalTraceDerivFactors (Array< OneD, Array< OneD, NekDouble > > &factors, Array< OneD, Array< OneD, NekDouble > > &d0factors, Array< OneD, Array< OneD, NekDouble > > &d1factors)
 
IndexMapValuesSharedPtr GetIndexMap (const IndexMapKey &ikey)
 
void AlignVectorToCollapsedDir (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
ExpansionSharedPtr GetLeftAdjacentElementExp () const
 
ExpansionSharedPtr GetRightAdjacentElementExp () const
 
int GetLeftAdjacentElementTrace () const
 
int GetRightAdjacentElementTrace () const
 
void SetAdjacentElementExp (int traceid, ExpansionSharedPtr &e)
 
StdRegions::Orientation GetTraceOrient (int trace)
 
void SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void DivideByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Divided by the metric jacobi and quadrature weights. More...
 
void GetTraceQFactors (const int trace, Array< OneD, NekDouble > &outarray)
 Extract the metric factors to compute the contravariant fluxes along edge edge and stores them into outarray following the local edge orientation (i.e. anticlockwise convention). More...
 
void GetTracePhysVals (const int trace, const StdRegions::StdExpansionSharedPtr &TraceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=StdRegions::eNoOrientation)
 
void GetTracePhysMap (const int edge, Array< OneD, int > &outarray)
 
void ReOrientTracePhysMap (const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1)
 
const NormalVectorGetTraceNormal (const int id)
 
void ComputeTraceNormal (const int id)
 
const Array< OneD, const NekDouble > & GetPhysNormals (void)
 
void SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
void SetUpPhysNormals (const int trace)
 
void AddRobinMassMatrix (const int traceid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
void TraceNormLen (const int traceid, NekDouble &h, NekDouble &p)
 
void AddRobinTraceContribution (const int traceid, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs)
 
const Array< OneD, const NekDouble > & GetElmtBndNormDirElmtLen (const int nbnd) const
 
void StdDerivBaseOnTraceMat (Array< OneD, DNekMatSharedPtr > &DerivMat)
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion
 StdExpansion ()
 Default Constructor. More...
 
 StdExpansion (const int numcoeffs, const int numbases, const LibUtilities::BasisKey &Ba=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bb=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bc=LibUtilities::NullBasisKey)
 Constructor. More...
 
 StdExpansion (const StdExpansion &T)
 Copy Constructor. More...
 
virtual ~StdExpansion ()
 Destructor. More...
 
int GetNumBases () const
 This function returns the number of 1D bases used in the expansion. More...
 
const Array< OneD, const LibUtilities::BasisSharedPtr > & GetBase () const
 This function gets the shared point to basis. More...
 
const LibUtilities::BasisSharedPtrGetBasis (int dir) const
 This function gets the shared point to basis in the dir direction. More...
 
int GetNcoeffs (void) const
 This function returns the total number of coefficients used in the expansion. More...
 
int GetTotPoints () const
 This function returns the total number of quadrature points used in the element. More...
 
LibUtilities::BasisType GetBasisType (const int dir) const
 This function returns the type of basis used in the dir direction. More...
 
int GetBasisNumModes (const int dir) const
 This function returns the number of expansion modes in the dir direction. More...
 
int EvalBasisNumModesMax (void) const
 This function returns the maximum number of expansion modes over all local directions. More...
 
LibUtilities::PointsType GetPointsType (const int dir) const
 This function returns the type of quadrature points used in the dir direction. More...
 
int GetNumPoints (const int dir) const
 This function returns the number of quadrature points in the dir direction. More...
 
const Array< OneD, const NekDouble > & GetPoints (const int dir) const
 This function returns a pointer to the array containing the quadrature points in dir direction. More...
 
int GetNverts () const
 This function returns the number of vertices of the expansion domain. More...
 
int GetTraceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th trace. More...
 
int GetTraceIntNcoeffs (const int i) const
 
int GetTraceNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th trace. More...
 
const LibUtilities::BasisKey GetTraceBasisKey (const int i, int k=-1) const
 This function returns the basis key belonging to the i-th trace. More...
 
LibUtilities::PointsKey GetTracePointsKey (const int i, int k=-1) const
 This function returns the basis key belonging to the i-th trace. More...
 
int NumBndryCoeffs (void) const
 
int NumDGBndryCoeffs (void) const
 
const LibUtilities::PointsKey GetNodalPointsKey () const
 This function returns the type of expansion Nodal point type if defined. More...
 
int GetNtraces () const
 Returns the number of trace elements connected to this element. More...
 
LibUtilities::ShapeType DetShapeType () const
 This function returns the shape of the expansion domain. More...
 
std::shared_ptr< StdExpansionGetStdExp () const
 
std::shared_ptr< StdExpansionGetLinStdExp (void) const
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion () const
 
bool IsNodalNonTensorialExp ()
 
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space. More...
 
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space. More...
 
void FwdTransBndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble Integral (const Array< OneD, const NekDouble > &inarray)
 This function integrates the specified function over the domain. More...
 
void FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 This function fills the array outarray with the mode-th mode of the expansion. More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 this function calculates the inner product of a given function f with the different modes of the expansion More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check)
 
void IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
int GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void SetElmtId (const int id)
 Set the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2=NullNekDouble1DArray, Array< OneD, NekDouble > &coords_3=NullNekDouble1DArray)
 this function returns the physical coordinates of the quadrature points of the expansion More...
 
void GetCoord (const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
 given the coordinates of a point of the element in the local collapsed coordinate system, this function calculates the physical coordinates of the point More...
 
DNekMatSharedPtr GetStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr GetStdStaticCondMatrix (const StdMatrixKey &mkey)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
void DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
int CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
NekDouble StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
int GetCoordim ()
 
void GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
void GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
int GetVertexMap (const int localVertexId, bool useCoeffPacking=false)
 
void GetTraceToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
 
void GetTraceCoeffMap (const unsigned int traceid, Array< OneD, unsigned int > &maparray)
 
void GetElmtTraceToTraceMap (const unsigned int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
 
void GetTraceInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eForwards)
 
void GetTraceNumModes (const int tid, int &numModes0, int &numModes1, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2)
 
void MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
DNekMatSharedPtr CreateGeneralMatrix (const StdMatrixKey &mkey)
 this function generates the mass matrix \(\mathbf{M}[i][j] = \int \phi_i(\mathbf{x}) \phi_j(\mathbf{x}) d\mathbf{x}\) More...
 
void GeneralMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
void ExponentialFilter (Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
 
void LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void 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)
 
void StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
 This function evaluates the first derivative of the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs, std::array< NekDouble, 6 > &secondOrderDerivs)
 
NekDouble PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluateBasis (const Array< OneD, const NekDouble > &coords, int mode)
 This function evaluates the basis function mode mode at a point coords of the domain. More...
 
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta. More...
 
void LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi)
 Convert local collapsed coordinates eta into local cartesian coordinate xi. More...
 
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual void v_DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
NekDouble Linf (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( L_\infty\) error \( |\epsilon|_\infty = \max |u - u_{exact}|\) where \( u_{exact}\) is given by the array sol. More...
 
NekDouble L2 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( L_2\) error, \( | \epsilon |_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 dx \right]^{1/2} d\xi_1 \) where \( u_{exact}\) is given by the array sol. More...
 
NekDouble H1 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( H^1\) error, \( | \epsilon |^1_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 + \nabla(u - u_{exact})\cdot\nabla(u - u_{exact})\cdot dx \right]^{1/2} d\xi_1 \) where \( u_{exact}\) is given by the array sol. More...
 
const LibUtilities::PointsKeyVector GetPointsKeys () const
 
DNekMatSharedPtr BuildInverseTransformationMatrix (const DNekScalMatSharedPtr &m_transformationmatrix)
 
void PhysInterpToSimplexEquiSpaced (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int npset=-1)
 This function performs an interpolation from the physical space points provided at input into an array of equispaced points which are not the collapsed coordinate. So for a tetrahedron you will only get a tetrahedral number of values. More...
 
void GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true)
 This function provides the connectivity of local simplices (triangles or tets) to connect the equispaced data points provided by PhysInterpToSimplexEquiSpaced. More...
 
void EquiSpacedToCoeffs (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs a projection/interpolation from the equispaced points sometimes used in post-processing onto the coefficient space. More...
 
template<class T >
std::shared_ptr< T > as ()
 
void IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
void GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 

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. More...
 
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase_SumFac (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GeneralMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
 
void LaplacianMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp_MatFree (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void 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 void v_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
 
virtual NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
virtual void v_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. More...
 
template<int DIR>
NekDouble BaryEvaluateBasis (const NekDouble &coord, const int &mode)
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals)
 Helper function to pass an unused value by reference into BaryEvaluate. More...
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals, NekDouble &deriv)
 

Protected Attributes

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

Detailed Description

Definition at line 71 of file Expansion.h.

Constructor & Destructor Documentation

◆ Expansion() [1/2]

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

Definition at line 43 of file Expansion.cpp.

45 std::bind(&Expansion::CreateIndexMap, this, std::placeholders::_1),
46 std::string("ExpansionIndexMap")),
47 m_geom(pGeom), m_metricinfo(m_geom->GetGeomFactors()),
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...
Definition: ErrorUtil.hpp:202
IndexMapValuesSharedPtr CreateIndexMap(const IndexMapKey &ikey)
Definition: Expansion.cpp:181
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:273
std::map< int, ExpansionWeakPtr > m_traceExp
Definition: Expansion.h:272
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:274
LibUtilities::NekManager< IndexMapKey, IndexMapValues, IndexMapKey::opLess > m_indexMapManager
Definition: Expansion.h:270
Array< OneD, LibUtilities::BasisSharedPtr > m_base
@ eDeformed
Geometry is curved or has non-constant factors.

References Nektar::SpatialDomains::eDeformed, Nektar::ErrorUtil::ewarning, 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)
Definition: Expansion.cpp:756

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)
Definition: Expansion.cpp:775

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)
Definition: Expansion.cpp:890

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)
Definition: Expansion.cpp:898

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)
Definition: Expansion.cpp:923

References v_AlignVectorToCollapsedDir().

Referenced by Nektar::LocalRegions::NodalTriExp::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)
Definition: Expansion.cpp:731

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)
Definition: Expansion.cpp:739

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 603 of file Expansion.cpp.

606{
607 int shapedim = dfdir.size();
608 int coordim = GetCoordim();
609 int nqtot = direction.size() / coordim;
610
611 for (int j = 0; j < shapedim; j++)
612 {
613 dfdir[j] = Array<OneD, NekDouble>(nqtot, 0.0);
614 for (int k = 0; k < coordim; k++)
615 {
616 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
617 {
618 Vmath::Vvtvp(nqtot, &df[shapedim * k + j][0], 1,
619 &direction[k * nqtot], 1, &dfdir[j][0], 1,
620 &dfdir[j][0], 1);
621 }
622 else
623 {
624 Vmath::Svtvp(nqtot, df[shapedim * k + j][0],
625 &direction[k * nqtot], 1, &dfdir[j][0], 1,
626 &dfdir[j][0], 1);
627 }
628 }
629 }
630}
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 454 of file Expansion.cpp.

455{
456 unsigned int nqtot = GetTotPoints();
457 SpatialDomains::GeomType type = m_metricinfo->GetGtype();
459 if (type == SpatialDomains::eRegular ||
461 {
463 Array<OneD, NekDouble>(nqtot, m_metricinfo->GetJac(p)[0]);
464 }
465 else
466 {
468 }
469
472}
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:134
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.

References Nektar::LocalRegions::eMetricQuadrature, Nektar::SpatialDomains::eMovingRegular, Nektar::SpatialDomains::eRegular, Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::GetTotPoints(), m_metricinfo, m_metrics, CellMLToNektar.cellml_metadata::p, 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)
Definition: Expansion.cpp:868

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)
Definition: ErrorUtil.hpp:208
#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)
Definition: StdExpansion.h:684
std::shared_ptr< IndexMapValues > IndexMapValuesSharedPtr
Definition: IndexMapKey.h:126

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 }
410 return returnval;
411}
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:265
DNekScalMatSharedPtr GetLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:84
void GetBoundaryMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:669
DNekBlkMatSharedPtr GetStdStaticCondMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:608
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:674
@ eNoGeomType
No type defined.
NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag > DNekScalMat
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:77
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:79
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble

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, and Nektar::StdRegions::StdExpansion::NumBndryCoeffs().

◆ 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)
Definition: Expansion.cpp:784

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)
Definition: Expansion.cpp:434

References v_DivideByQuadratureMetric().

◆ DropLocMatrix()

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

Definition at line 90 of file Expansion.cpp.

91{
92 return v_DropLocMatrix(mkey);
93}
virtual void v_DropLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:413

References v_DropLocMatrix().

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

◆ 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)
Definition: Expansion.cpp:746

References v_ExtractDataToCoeffs().

◆ GetElmtBndNormDirElmtLen()

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

Definition at line 914 of file Expansion.cpp.

916{
917 auto x = m_elmtBndNormDirElmtLen.find(nbnd);
919 "m_elmtBndNormDirElmtLen normal not computed.");
920 return x->second;
921}
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:286

References ASSERTL0, and m_elmtBndNormDirElmtLen.

◆ GetGeom()

SpatialDomains::GeometrySharedPtr 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 441 of file Expansion.h.

442{
443 ASSERTL1(m_elementLeft.lock().get(), "Left adjacent element not set.");
444 return m_elementLeft.lock();
445}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
ExpansionWeakPtr m_elementLeft
Definition: Expansion.h:277

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.
Definition: StdExpansion.h:367
StdRegions::ConstFactorMap factors

References Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::VarcoeffHashingTest::factors, 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 633 of file Expansion.cpp.

635{
636 int coordim = GetCoordim();
637
638 int nquad0, nquad1, nquad2;
639 int nqtot = 1;
640 switch (shapedim)
641 {
642 case 2:
643 {
644 nquad0 = m_base[0]->GetNumPoints();
645 nquad1 = m_base[1]->GetNumPoints();
646 nqtot = nquad0 * nquad1;
647 break;
648 }
649 case 3:
650 {
651 nquad0 = m_base[0]->GetNumPoints();
652 nquad1 = m_base[1]->GetNumPoints();
653 nquad2 = m_base[2]->GetNumPoints();
654 nqtot = nquad0 * nquad1 * nquad2;
655 break;
656 }
657 default:
658 break;
659 }
660
661 StdRegions::VarCoeffType MMFCoeffs[15] = {
670
671 Array<OneD, NekDouble> MF(coordim * nqtot);
672 Array<OneD, NekDouble> tmp(nqtot);
673 for (int k = 0; k < coordim; k++)
674 {
675 StdRegions::VarCoeffMap::const_iterator MFdir =
676 varcoeffs.find(MMFCoeffs[5 * dir + k]);
677 tmp = MFdir->second.GetValue();
678
679 Vmath::Vcopy(nqtot, &tmp[0], 1, &MF[k * nqtot], 1);
680 }
681
682 return MF;
683}
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 686 of file Expansion.cpp.

688{
689 int indxDiv = 3;
690
691 StdRegions::VarCoeffType MMFCoeffs[15] = {
700
701 StdRegions::VarCoeffMap::const_iterator MFdir =
702 varcoeffs.find(MMFCoeffs[5 * dir + indxDiv]);
703 Array<OneD, NekDouble> MFDiv = MFdir->second.GetValue();
704
705 return MFDiv;
706}

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 709 of file Expansion.cpp.

711{
712 int indxMag = 4;
713
714 StdRegions::VarCoeffType MMFCoeffs[15] = {
723
724 StdRegions::VarCoeffMap::const_iterator MFdir =
725 varcoeffs.find(MMFCoeffs[5 * dir + indxMag]);
726 Array<OneD, NekDouble> MFmag = MFdir->second.GetValue();
727
728 return MFmag;
729}

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()
Definition: Expansion.cpp:873

References v_GetPhysNormals().

◆ GetRightAdjacentElementExp()

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

Definition at line 447 of file Expansion.h.

448{
449 ASSERTL1(m_elementLeft.lock().get(), "Right adjacent element not set.");
450
451 return m_elementRight.lock();
452}
ExpansionWeakPtr m_elementRight
Definition: Expansion.h:278

References ASSERTL1, m_elementLeft, and m_elementRight.

◆ GetRightAdjacentElementTrace()

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

Definition at line 459 of file Expansion.h.

460{
461 return m_elementTraceRight;
462}

References m_elementTraceRight.

◆ GetTraceExp()

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

Definition at line 414 of file Expansion.h.

415{
416 ASSERTL1(traceid < GetNtraces(), "Trace is out of range.");
417
418 ExpansionSharedPtr returnval;
419
420 if (m_traceExp.count(traceid))
421 {
422 // Use stored value
423 returnval = m_traceExp[traceid].lock();
424 }
425 else
426 {
427 // Generate trace exp
428 v_GenTraceExp(traceid, returnval);
429 }
430
431 return returnval;
432}
virtual void v_GenTraceExp(const int traceid, ExpansionSharedPtr &exp)
Definition: Expansion.cpp:861
int GetNtraces() const
Returns the number of trace elements connected to this element.
Definition: StdExpansion.h:351
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)
Definition: Expansion.cpp:845

References v_GetTracePhysMap().

Referenced by 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)
Definition: Expansion.cpp:834

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)
Definition: Expansion.cpp:826

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)
Definition: Expansion.cpp:803

References Nektar::VarcoeffHashingTest::factors, and v_NormalTraceDerivFactors().

◆ 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)
Definition: Expansion.cpp:852

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
178 m_metricinfo = m_geom->GetGeomFactors();
179}

References m_geom, m_metricinfo, and m_metrics.

◆ SetAdjacentElementExp()

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

Definition at line 464 of file Expansion.h.

466{
467 if (m_elementLeft.lock().get())
468 {
469 m_elementRight = exp;
470 m_elementTraceRight = traceid;
471 }
472 else
473 {
474 m_elementLeft = exp;
475 m_elementTraceLeft = traceid;
476 }
477}

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
Definition: Expansion.cpp:818

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)
Definition: Expansion.cpp:879

References v_SetPhysNormals().

◆ SetTraceExp()

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

Definition at line 434 of file Expansion.h.

435{
436 ASSERTL1(traceid < GetNtraces(), "Trace out of range.");
437
438 m_traceExp[traceid] = exp;
439}

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)
Definition: Expansion.cpp:885

References v_SetUpPhysNormals().

◆ StdDerivBaseOnTraceMat()

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

Definition at line 474 of file Expansion.cpp.

475{
476 int nquad = GetTotPoints();
477 int ntraces = GetNtraces();
478 int ndir = m_base.size();
479
480 Array<OneD, NekDouble> coeffs(m_ncoeffs);
481 Array<OneD, NekDouble> phys(nquad);
482
483 Array<OneD, Array<OneD, int>> traceids(ntraces);
484
485 int tottracepts = 0;
486 for (int i = 0; i < ntraces; ++i)
487 {
488 GetTracePhysMap(i, traceids[i]);
489 tottracepts += GetTraceNumPoints(i);
490 }
491
492 // initialise array to null so can call for
493 // differnt dimensions
494 Array<OneD, Array<OneD, NekDouble>> Deriv(3, NullNekDouble1DArray);
495
496 DerivMat = Array<OneD, DNekMatSharedPtr>(ndir);
497
498 for (int i = 0; i < ndir; ++i)
499 {
500 Deriv[i] = Array<OneD, NekDouble>(nquad);
501 DerivMat[i] =
503 }
504
505 for (int i = 0; i < m_ncoeffs; ++i)
506 {
507 Vmath::Zero(m_ncoeffs, coeffs, 1);
508 coeffs[i] = 1.0;
509 BwdTrans(coeffs, phys);
510
511 // dphi_i/d\xi_1, dphi_i/d\xi_2 dphi_i/d\xi_3
512 StdPhysDeriv(phys, Deriv[0], Deriv[1], Deriv[2]);
513
514 int cnt = 0;
515 for (int j = 0; j < ntraces; ++j)
516 {
517 int nTracePts = GetTraceNumPoints(j);
518 for (int k = 0; k < nTracePts; ++k)
519 {
520 for (int d = 0; d < ndir; ++d)
521 {
522 (*DerivMat[d])(i, cnt + k) = Deriv[d][traceids[j][k]];
523 }
524 }
525 cnt += nTracePts;
526 }
527 }
528}
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.
Definition: StdExpansion.h:281
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Backward transformation from coefficient space to physical space.
Definition: StdExpansion.h:424
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)
Definition: StdExpansion.h:882
std::vector< double > d(NPUPPER *NPUPPER)
static Array< OneD, NekDouble > NullNekDouble1DArray
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::UnitTests::d(), 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)
Definition: Expansion.cpp:907

References CellMLToNektar.cellml_metadata::p, and 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 766 of file Expansion.cpp.

771{
772 NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
773}

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 756 of file Expansion.cpp.

762{
763 NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
764}

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

Referenced by 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 775 of file Expansion.cpp.

780{
781 NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
782}

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 890 of file Expansion.cpp.

894{
895 NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
896}

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 898 of file Expansion.cpp.

903{
904 NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
905}

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 731 of file Expansion.cpp.

734{
735 NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
737}
static DNekMatSharedPtr NullDNekMatSharedPtr
Definition: NekTypeDefs.hpp:83

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 739 of file Expansion.cpp.

741{
742 NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
744}

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 784 of file Expansion.cpp.

790{
791 NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
792}

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 434 of file Expansion.cpp.

437{
438 const int nqtot = GetTotPoints();
439
440 if (m_metrics.count(eMetricQuadrature) == 0)
441 {
443 }
444
445 Vmath::Vdiv(nqtot, inarray, 1, m_metrics[eMetricQuadrature], 1, outarray,
446 1);
447}
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 861 of file Expansion.cpp.

863{
865 "Method does not exist for this shape or library");
866}

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 315 of file Expansion.h.

316 {
317 return m_geom->GetCoordim();
318 }

References 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 530 of file Expansion.cpp.

533{
534 ASSERTL1(m_geom, "m_geom not defined");
535
536 // get physical points defined in Geom
537 m_geom->FillGeom();
538
539 const int expDim = m_base.size();
540 int nqGeom = 1;
541 bool doCopy = true;
542
543 Array<OneD, LibUtilities::BasisSharedPtr> CBasis(expDim);
544 Array<OneD, Array<OneD, NekDouble>> tmp(3);
545
546 for (int i = 0; i < expDim; ++i)
547 {
548 CBasis[i] = m_geom->GetXmap()->GetBasis(i);
549 nqGeom *= CBasis[i]->GetNumPoints();
550 doCopy = doCopy &&
551 m_base[i]->GetBasisKey().SamePoints(CBasis[i]->GetBasisKey());
552 }
553
554 tmp[0] = coords_0;
555 tmp[1] = coords_1;
556 tmp[2] = coords_2;
557
558 if (doCopy)
559 {
560 for (int i = 0; i < m_geom->GetCoordim(); ++i)
561 {
562 m_geom->GetXmap()->BwdTrans(m_geom->GetCoeffs(i), tmp[i]);
563 }
564 }
565 else
566 {
567 for (int i = 0; i < m_geom->GetCoordim(); ++i)
568 {
569 Array<OneD, NekDouble> tmpGeom(nqGeom);
570 m_geom->GetXmap()->BwdTrans(m_geom->GetCoeffs(i), tmpGeom);
571
572 switch (expDim)
573 {
574 case 1:
575 {
577 CBasis[0]->GetPointsKey(), &tmpGeom[0],
578 m_base[0]->GetPointsKey(), &tmp[i][0]);
579 break;
580 }
581 case 2:
582 {
584 CBasis[0]->GetPointsKey(), CBasis[1]->GetPointsKey(),
585 &tmpGeom[0], m_base[0]->GetPointsKey(),
586 m_base[1]->GetPointsKey(), &tmp[i][0]);
587 break;
588 }
589 case 3:
590 {
592 CBasis[0]->GetPointsKey(), CBasis[1]->GetPointsKey(),
593 CBasis[2]->GetPointsKey(), &tmpGeom[0],
594 m_base[0]->GetPointsKey(), m_base[1]->GetPointsKey(),
595 m_base[2]->GetPointsKey(), &tmp[i][0]);
596 break;
597 }
598 }
599 }
600 }
601}
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::LibUtilities::Interp1D(), Nektar::LibUtilities::Interp2D(), Nektar::LibUtilities::Interp3D(), Nektar::StdRegions::StdExpansion::m_base, and m_geom.

Referenced by Nektar::LocalRegions::NodalTriExp::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 873 of file Expansion.cpp.

874{
875 NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
877}

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

Definition at line 834 of file Expansion.cpp.

840{
842 "Method does not exist for this shape or library");
843}

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 826 of file Expansion.cpp.

829{
831 "Method does not exist for this shape or library");
832}

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 419 of file Expansion.cpp.

422{
423 const int nqtot = GetTotPoints();
424
425 if (m_metrics.count(eMetricQuadrature) == 0)
426 {
428 }
429
430 Vmath::Vmul(nqtot, m_metrics[eMetricQuadrature], 1, inarray, 1, outarray,
431 1);
432}
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 803 of file Expansion.cpp.

807{
809 "This function is only valid for "
810 "shape expansion in LocalRegions, not parant class");
811}

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 852 of file Expansion.cpp.

856{
858 "Method does not exist for this shape or library");
859}

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 818 of file Expansion.cpp.

822{
823 NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
824}

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 879 of file Expansion.cpp.

881{
882 NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
883}

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 885 of file Expansion.cpp.

886{
887 NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
888}

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 907 of file Expansion.cpp.

910{
911 NEKERROR(ErrorUtil::efatal, "This method has not been defined");
912}

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 794 of file Expansion.cpp.

796{
798 "This function is only valid for "
799 "shape expansion in LocalRegions, not parant class");
800 return 0.0;
801}

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)
Definition: Expansion.cpp:794

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 278 of file Expansion.h.

Referenced by GetRightAdjacentElementExp(), and SetAdjacentElementExp().

◆ m_elementTraceLeft

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

Definition at line 279 of file Expansion.h.

Referenced by GetLeftAdjacentElementTrace(), and SetAdjacentElementExp().

◆ m_elementTraceRight

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

Definition at line 280 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::GeometrySharedPtr Nektar::LocalRegions::Expansion::m_geom
protected

Definition at line 273 of file Expansion.h.

Referenced by Nektar::LocalRegions::SegExp::CreateMatrix(), Expansion(), Nektar::LocalRegions::NodalTriExp::GetCoord(), Nektar::LocalRegions::PointExp::GetCoords(), GetGeom(), Nektar::LocalRegions::PointExp::GetGeom(), Nektar::LocalRegions::Expansion0D::GetGeom0D(), Nektar::LocalRegions::Expansion1D::GetGeom1D(), Nektar::LocalRegions::Expansion2D::GetGeom2D(), Nektar::LocalRegions::Expansion3D::GetGeom3D(), Nektar::LocalRegions::NodalTriExp::PhysEvaluate(), 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::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::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(), Nektar::LocalRegions::TetExp::v_PhysEvaluate(), and Nektar::LocalRegions::NodalTriExp::v_PhysEvaluate().

◆ m_indexMapManager

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

Definition at line 270 of file Expansion.h.

Referenced by GetIndexMap().

◆ m_metricinfo

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

Definition at line 274 of file Expansion.h.

Referenced by ComputeGmatcdotMF(), ComputeQuadratureMetric(), Nektar::LocalRegions::Expansion2D::CreateMatrix(), Nektar::LocalRegions::Expansion3D::CreateMatrix(), Nektar::LocalRegions::SegExp::CreateMatrix(), CreateStaticCondMatrix(), Expansion(), GetMetricInfo(), Nektar::LocalRegions::NodalTriExp::Integral(), Nektar::LocalRegions::NodalTriExp::PhysDeriv(), 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::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::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