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)
 
virtual ~Expansion ()
 
void SetTraceExp (const int traceid, ExpansionSharedPtr &f)
 
ExpansionSharedPtr GetTraceExp (const int traceid)
 
DNekScalMatSharedPtr GetLocMatrix (const LocalRegions::MatrixKey &mkey)
 
DNekScalMatSharedPtr GetLocMatrix (const StdRegions::MatrixType mtype, const StdRegions::ConstFactorMap &factors=StdRegions::NullConstFactorMap, const StdRegions::VarCoeffMap &varcoeffs=StdRegions::NullVarCoeffMap)
 
SpatialDomains::GeometrySharedPtr GetGeom () const
 
void Reset ()
 
IndexMapValuesSharedPtr CreateIndexMap (const IndexMapKey &ikey)
 
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)
 
virtual void AddRobinTraceContribution (const int traceid, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs)
 
const Array< OneD, const NekDouble > & GetElmtBndNormDirElmtLen (const int nbnd) const
 
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 (void) const
 
std::shared_ptr< StdExpansionGetLinStdExp (void) const
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion ()
 
bool IsNodalNonTensorialExp ()
 
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space. More...
 
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space. More...
 
void FwdTransBndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble Integral (const Array< OneD, const NekDouble > &inarray)
 This function integrates the specified function over the domain. More...
 
void FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 This function fills the array outarray with the mode-th mode of the expansion. More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 this function calculates the inner product of a given function f with the different modes of the expansion More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check)
 
void IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
int GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void SetElmtId (const int id)
 Set the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2=NullNekDouble1DArray, Array< OneD, NekDouble > &coords_3=NullNekDouble1DArray)
 this function returns the physical coordinates of the quadrature points of the expansion More...
 
void GetCoord (const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
 given the coordinates of a point of the element in the local collapsed coordinate system, this function calculates the physical coordinates of the point More...
 
DNekMatSharedPtr GetStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr GetStdStaticCondMatrix (const StdMatrixKey &mkey)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble >> &Fvec, Array< OneD, NekDouble > &outarray)
 
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
void DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
int CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
NekDouble StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
int GetCoordim ()
 
void GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
void GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
int GetVertexMap (const int localVertexId, bool useCoeffPacking=false)
 
void GetTraceToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
 
void GetTraceCoeffMap (const unsigned int traceid, Array< OneD, unsigned int > &maparray)
 
void GetElmtTraceToTraceMap (const unsigned int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
 
void GetTraceInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eForwards)
 
void GetTraceNumModes (const int tid, int &numModes0, int &numModes1, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2)
 
void MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
DNekMatSharedPtr CreateGeneralMatrix (const StdMatrixKey &mkey)
 this function generates the mass matrix \(\mathbf{M}[i][j] = \int \phi_i(\mathbf{x}) \phi_j(\mathbf{x}) d\mathbf{x}\) More...
 
void GeneralMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
void ExponentialFilter (Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
 
void LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
DNekMatSharedPtr GenMatrix (const StdMatrixKey &mkey)
 
void PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void PhysDeriv_s (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds)
 
void PhysDeriv_n (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn)
 
void PhysDirectionalDeriv (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &outarray)
 
void StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluate (const Array< OneD, 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)
 
virtual void v_MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_DivideByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_ComputeLaplacianMetric ()
 
virtual void v_GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3)
 
Array< OneD, NekDoublev_GetMF (const int dir, const int shapedim, const StdRegions::VarCoeffMap &varcoeffs)
 
Array< OneD, NekDoublev_GetMFDiv (const int dir, const StdRegions::VarCoeffMap &varcoeffs)
 
Array< OneD, NekDoublev_GetMFMag (const int dir, const StdRegions::VarCoeffMap &varcoeffs)
 
virtual DNekScalMatSharedPtr v_GetLocMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual 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)
 
virtual void v_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_GetTraceQFactors (const int trace, Array< OneD, NekDouble > &outarray)
 
virtual void v_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 (void)
 
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_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase_SumFac (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GeneralMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
 
void LaplacianMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp_MatFree (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void HelmholtzMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
 
virtual NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
virtual void v_GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 
virtual void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
template<int DIR>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals)
 This function performs the barycentric interpolation of the polynomial stored in coord at a point physvals using barycentric interpolation weights in direction. More...
 
template<int DIR>
NekDouble BaryEvaluateBasis (const NekDouble &coord, const int &mode)
 

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

Constructor & Destructor Documentation

◆ Expansion() [1/2]

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

Definition at line 47 of file Expansion.cpp.

49  std::bind(&Expansion::CreateIndexMap, this, std::placeholders::_1),
50  std::string("ExpansionIndexMap")),
51  m_geom(pGeom), m_metricinfo(m_geom->GetGeomFactors()),
53 {
54  if (!m_metricinfo)
55  {
56  return;
57  }
58 
59  if (!m_metricinfo->IsValid())
60  {
61  int nDim = m_base.size();
62  string type = "regular";
63  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
64  {
65  type = "deformed";
66  }
67 
68  stringstream err;
69  err << nDim << "D " << type << " Jacobian not positive "
70  << "(element ID = " << m_geom->GetGlobalID() << ") "
71  << "(first vertex ID = " << m_geom->GetVid(0) << ")";
72  NEKERROR(ErrorUtil::ewarning, err.str());
73  }
74 
75  m_traceExp.empty();
76 }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
IndexMapValuesSharedPtr CreateIndexMap(const IndexMapKey &ikey)
Definition: Expansion.cpp:180
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:272
std::map< int, ExpansionWeakPtr > m_traceExp
Definition: Expansion.h:271
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:273
LibUtilities::NekManager< IndexMapKey, IndexMapValues, IndexMapKey::opLess > m_indexMapManager
Definition: Expansion.h:269
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 78 of file Expansion.cpp.

79  : StdExpansion(pSrc), m_indexMapManager(pSrc.m_indexMapManager),
80  m_geom(pSrc.m_geom), m_metricinfo(pSrc.m_metricinfo)
81 {
82 }
StdExpansion()
Default Constructor.

◆ ~Expansion()

Nektar::LocalRegions::Expansion::~Expansion ( )
virtual

Definition at line 84 of file Expansion.cpp.

85 {
86 }

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

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

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

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

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

131 {
132  v_AddFaceNormBoundaryInt(face, FaceExp, Fn, outarray);
133 }
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:767

References v_AddFaceNormBoundaryInt().

◆ AddRobinMassMatrix()

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

Definition at line 242 of file Expansion.h.

245  {
246  v_AddRobinMassMatrix(traceid, primCoeffs, inoutmat);
247  }
virtual void v_AddRobinMassMatrix(const int face, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
Definition: Expansion.cpp:885

References v_AddRobinMassMatrix().

◆ AddRobinTraceContribution()

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

Definition at line 254 of file Expansion.h.

257  {
258  v_AddRobinTraceContribution(traceid, primCoeffs, incoeffs, coeffs);
259  }
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:893

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

153  {
154  v_AlignVectorToCollapsedDir(dir, inarray, outarray);
155  }
virtual void v_AlignVectorToCollapsedDir(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
Definition: Expansion.cpp:916

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

96 {
97  return v_BuildTransformationMatrix(r_bnd, matrixType);
98 }
virtual DNekMatSharedPtr v_BuildTransformationMatrix(const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
Definition: Expansion.cpp:725

References v_BuildTransformationMatrix().

◆ BuildVertexMatrix()

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

Definition at line 100 of file Expansion.cpp.

101 {
102  return v_BuildVertexMatrix(r_bnd);
103 }
virtual DNekMatSharedPtr v_BuildVertexMatrix(const DNekScalMatSharedPtr &r_bnd)
Definition: Expansion.cpp:733

References v_BuildVertexMatrix().

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

600 {
601  int shapedim = dfdir.size();
602  int coordim = GetCoordim();
603  int nqtot = direction.size() / coordim;
604 
605  for (int j = 0; j < shapedim; j++)
606  {
607  dfdir[j] = Array<OneD, NekDouble>(nqtot, 0.0);
608  for (int k = 0; k < coordim; k++)
609  {
610  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
611  {
612  Vmath::Vvtvp(nqtot, &df[shapedim * k + j][0], 1,
613  &direction[k * nqtot], 1, &dfdir[j][0], 1,
614  &dfdir[j][0], 1);
615  }
616  else
617  {
618  Vmath::Svtvp(nqtot, df[shapedim * k + j][0],
619  &direction[k * nqtot], 1, &dfdir[j][0], 1,
620  &dfdir[j][0], 1);
621  }
622  }
623  }
624 }
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
Definition: Vmath.cpp:622
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:574

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

Referenced by Nektar::LocalRegions::Expansion2D::CreateMatrix(), Nektar::LocalRegions::HexExp::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 448 of file Expansion.cpp.

449 {
450  unsigned int nqtot = GetTotPoints();
451  SpatialDomains::GeomType type = m_metricinfo->GetGtype();
453  if (type == SpatialDomains::eRegular ||
455  {
457  Array<OneD, NekDouble>(nqtot, m_metricinfo->GetJac(p)[0]);
458  }
459  else
460  {
462  }
463 
466 }
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140
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:250
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 222 of file Expansion.h.

223  {
225  }
virtual void v_ComputeTraceNormal(const int id)
Definition: Expansion.cpp:861

References v_ComputeTraceNormal().

◆ CreateIndexMap()

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

Definition at line 180 of file Expansion.cpp.

181 {
182  IndexMapValuesSharedPtr returnval;
183 
184  IndexMapType itype = ikey.GetIndexMapType();
185 
186  int entity = ikey.GetIndexEntity();
187 
188  StdRegions::Orientation orient = ikey.GetIndexOrientation();
189 
190  Array<OneD, unsigned int> map;
191  Array<OneD, int> sign;
192 
193  switch (itype)
194  {
195  case eEdgeToElement:
196  {
197  GetTraceToElementMap(entity, map, sign, orient);
198  }
199  break;
200  case eFaceToElement:
201  {
202  GetTraceToElementMap(entity, map, sign, orient);
203  }
204  break;
205  case eEdgeInterior:
206  {
207  ASSERTL0(false, "Boundary Index Map not implemented yet.");
208  // v_GetEdgeInteriorMap(entity,orient,map,sign);
209  }
210  break;
211  case eFaceInterior:
212  {
213  ASSERTL0(false, "Boundary Index Map not implemented yet.");
214  // v_GetFaceInteriorMap(entity,orient,map,sign);
215  }
216  break;
217  case eBoundary:
218  {
219  ASSERTL0(false, "Boundary Index Map not implemented yet.");
220  }
221  break;
222  case eVertex:
223  {
224  ASSERTL0(false, "Vertex Index Map not implemented yet.");
225  }
226  break;
227  default:
228  {
229  ASSERTL0(false, "The Index Map you are requiring "
230  "is not between the possible options.");
231  }
232  }
233 
234  returnval = MemoryManager<IndexMapValues>::AllocateSharedPtr(map.size());
235 
236  for (int i = 0; i < map.size(); i++)
237  {
238  (*returnval)[i].index = map[i];
239  (*returnval)[i].sign = sign[i];
240  }
241 
242  return returnval;
243 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:15
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:692
std::shared_ptr< IndexMapValues > IndexMapValuesSharedPtr
Definition: IndexMapKey.h:128

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
293  case StdRegions::eMass:
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);
345  DNekMatSharedPtr B =
347  DNekMatSharedPtr C =
349  DNekMatSharedPtr D =
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:272
DNekScalMatSharedPtr GetLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:88
void GetBoundaryMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:677
DNekBlkMatSharedPtr GetStdStaticCondMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:616
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:682
@ 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 135 of file Expansion.cpp.

140 {
141  v_DGDeriv(dir, inarray, EdgeExp, coeffs, outarray);
142 }
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:775

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

183  {
184  v_DivideByQuadratureMetric(inarray, outarray);
185  }
virtual void v_DivideByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:428

References v_DivideByQuadratureMetric().

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

109 {
110  v_ExtractDataToCoeffs(data, nummodes, nmodes_offset, coeffs, fromType);
111 }
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:741

References v_ExtractDataToCoeffs().

◆ GetElmtBndNormDirElmtLen()

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

Definition at line 907 of file Expansion.cpp.

909 {
910  auto x = m_elmtBndNormDirElmtLen.find(nbnd);
912  "m_elmtBndNormDirElmtLen normal not computed.");
913  return x->second;
914 }
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:285

References ASSERTL0, and m_elmtBndNormDirElmtLen.

◆ GetGeom()

SpatialDomains::GeometrySharedPtr Nektar::LocalRegions::Expansion::GetGeom ( ) const

◆ GetIndexMap()

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

Definition at line 145 of file Expansion.h.

146  {
147  return m_indexMapManager[ikey];
148  }

References m_indexMapManager.

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

◆ GetLeftAdjacentElementExp()

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

Definition at line 433 of file Expansion.h.

434 {
435  ASSERTL1(m_elementLeft.lock().get(), "Left adjacent element not set.");
436  return m_elementLeft.lock();
437 }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
ExpansionWeakPtr m_elementLeft
Definition: Expansion.h:276

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(), 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 157 of file Expansion.cpp.

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

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

◆ GetMetricInfo()

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

◆ GetPhysNormals()

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

Definition at line 227 of file Expansion.h.

228  {
229  return v_GetPhysNormals();
230  }
virtual const Array< OneD, const NekDouble > & v_GetPhysNormals(void)
Definition: Expansion.cpp:867

References v_GetPhysNormals().

◆ GetRightAdjacentElementExp()

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

Definition at line 439 of file Expansion.h.

440 {
441  ASSERTL1(m_elementLeft.lock().get(), "Right adjacent element not set.");
442 
443  return m_elementRight.lock();
444 }
ExpansionWeakPtr m_elementRight
Definition: Expansion.h:277

References ASSERTL1, m_elementLeft, and m_elementRight.

◆ GetRightAdjacentElementTrace()

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

Definition at line 451 of file Expansion.h.

452 {
453  return m_elementTraceRight;
454 }

References m_elementTraceRight.

◆ GetTraceExp()

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

Definition at line 406 of file Expansion.h.

407 {
408  ASSERTL1(traceid < GetNtraces(), "Trace is out of range.");
409 
410  ExpansionSharedPtr returnval;
411 
412  if (m_traceExp.count(traceid))
413  {
414  // Use stored value
415  returnval = m_traceExp[traceid].lock();
416  }
417  else
418  {
419  // Generate trace exp
420  v_GenTraceExp(traceid, returnval);
421  }
422 
423  return returnval;
424 }
virtual void v_GenTraceExp(const int traceid, ExpansionSharedPtr &exp)
Definition: Expansion.cpp:854
int GetNtraces() const
Returns the number of trace elements connected to this element.
Definition: StdExpansion.h:359
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68

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

Referenced by Nektar::LocalRegions::Expansion2D::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 208 of file Expansion.h.

209  {
210  v_GetTracePhysMap(edge, outarray);
211  }
virtual void v_GetTracePhysMap(const int edge, Array< OneD, int > &outarray)
Definition: Expansion.cpp:838

References v_GetTracePhysMap().

Referenced by StdDerivBaseOnTraceMat().

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

204  {
205  v_GetTracePhysVals(trace, TraceExp, inarray, outarray, orient);
206  }
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:828

References v_GetTracePhysVals().

Referenced by Nektar::LocalRegions::Expansion2D::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 193 of file Expansion.h.

195  {
196  v_GetTraceQFactors(trace, outarray);
197  }
virtual void v_GetTraceQFactors(const int trace, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:820

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

153 {
154  return v_NormalTraceDerivFactors(factors, d0factors, d1factors);
155 }
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:795

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

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

References v_ReOrientTracePhysMap().

◆ Reset()

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

Definition at line 171 of file Expansion.cpp.

172 {
173  // Clear metrics
174  m_metrics.clear();
175 
176  // Regenerate geometry factors
177  m_metricinfo = m_geom->GetGeomFactors();
178 }

References m_geom, m_metricinfo, and m_metrics.

◆ SetAdjacentElementExp()

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

Definition at line 456 of file Expansion.h.

458 {
459  if (m_elementLeft.lock().get())
460  {
461  m_elementRight = exp;
462  m_elementTraceRight = traceid;
463  }
464  else
465  {
466  m_elementLeft = exp;
467  m_elementTraceLeft = traceid;
468  }
469 }

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

175  {
176  v_SetCoeffsToOrientation(dir, inarray, outarray);
177  }
virtual void v_SetCoeffsToOrientation(StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: Expansion.cpp:812

References v_SetCoeffsToOrientation().

◆ SetPhysNormals()

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

Definition at line 232 of file Expansion.h.

233  {
234  v_SetPhysNormals(normal);
235  }
virtual void v_SetPhysNormals(Array< OneD, const NekDouble > &normal)
Definition: Expansion.cpp:873

References v_SetPhysNormals().

◆ SetTraceExp()

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

Definition at line 426 of file Expansion.h.

427 {
428  ASSERTL1(traceid < GetNtraces(), "Trace out of range.");
429 
430  m_traceExp[traceid] = exp;
431 }

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

◆ SetUpPhysNormals()

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

Definition at line 237 of file Expansion.h.

238  {
239  v_SetUpPhysNormals(trace);
240  }
virtual void v_SetUpPhysNormals(const int id)
Definition: Expansion.cpp:879

References v_SetUpPhysNormals().

◆ StdDerivBaseOnTraceMat()

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

Definition at line 468 of file Expansion.cpp.

469 {
470  int nquad = GetTotPoints();
471  int ntraces = GetNtraces();
472  int ndir = m_base.size();
473 
474  Array<OneD, NekDouble> coeffs(m_ncoeffs);
475  Array<OneD, NekDouble> phys(nquad);
476 
477  Array<OneD, Array<OneD, int>> traceids(ntraces);
478 
479  int tottracepts = 0;
480  for (int i = 0; i < ntraces; ++i)
481  {
482  GetTracePhysMap(i, traceids[i]);
483  tottracepts += GetTraceNumPoints(i);
484  }
485 
486  // initialise array to null so can call for
487  // differnt dimensions
488  Array<OneD, Array<OneD, NekDouble>> Deriv(3, NullNekDouble1DArray);
489 
490  DerivMat = Array<OneD, DNekMatSharedPtr>(ndir);
491 
492  for (int i = 0; i < ndir; ++i)
493  {
494  Deriv[i] = Array<OneD, NekDouble>(nquad);
495  DerivMat[i] =
497  }
498 
499  for (int i = 0; i < m_ncoeffs; ++i)
500  {
501  Vmath::Zero(m_ncoeffs, coeffs, 1);
502  coeffs[i] = 1.0;
503  BwdTrans(coeffs, phys);
504 
505  // dphi_i/d\xi_1, dphi_i/d\xi_2 dphi_i/d\xi_3
506  StdPhysDeriv(phys, Deriv[0], Deriv[1], Deriv[2]);
507 
508  int cnt = 0;
509  for (int j = 0; j < ntraces; ++j)
510  {
511  int nTracePts = GetTraceNumPoints(j);
512  for (int k = 0; k < nTracePts; ++k)
513  {
514  for (int d = 0; d < ndir; ++d)
515  {
516  (*DerivMat[d])(i, cnt + k) = Deriv[d][traceids[j][k]];
517  }
518  }
519  cnt += nTracePts;
520  }
521  }
522 }
void GetTracePhysMap(const int edge, Array< OneD, int > &outarray)
Definition: Expansion.h:208
int GetTraceNumPoints(const int i) const
This function returns the number of quadrature points belonging to the i-th trace.
Definition: StdExpansion.h:289
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:432
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:883
static Array< OneD, NekDouble > NullNekDouble1DArray
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492

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

◆ TraceNormLen()

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

Definition at line 249 of file Expansion.h.

250  {
251  v_TraceNormLen(traceid, h, p);
252  }
virtual void v_TraceNormLen(const int traceid, NekDouble &h, NekDouble &p)
Definition: Expansion.cpp:901

References CellMLToNektar.cellml_metadata::p, and v_TraceNormLen().

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

762 {
763  boost::ignore_unused(edge, EdgeExp, Fn, outarray);
764  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
765 }

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

754 {
755  boost::ignore_unused(edge, EdgeExp, Fx, Fy, outarray);
756  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
757 }

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

770 {
771  boost::ignore_unused(face, FaceExp, Fn, outarray);
772  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
773 }

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::Expansion1D, Nektar::LocalRegions::Expansion3D, and Nektar::LocalRegions::Expansion2D.

Definition at line 885 of file Expansion.cpp.

888 {
889  boost::ignore_unused(trace, primCoeffs, inoutmat);
890  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
891 }

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.

Definition at line 893 of file Expansion.cpp.

896 {
897  boost::ignore_unused(traceid, primCoeffs, incoeffs, coeffs);
898  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
899 }

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

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

Definition at line 916 of file Expansion.cpp.

919 {
920  boost::ignore_unused(dir, inarray, outarray);
921  NEKERROR(ErrorUtil::efatal, "v_AlignVectorToCollapsedDir is not defined");
922 }

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

Referenced by AlignVectorToCollapsedDir().

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

727 {
728  boost::ignore_unused(r_bnd, matrixType);
729  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
730  return NullDNekMatSharedPtr;
731 }
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::Expansion3D, and Nektar::LocalRegions::Expansion2D.

Definition at line 733 of file Expansion.cpp.

735 {
736  boost::ignore_unused(r_bnd);
737  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
738  return NullDNekMatSharedPtr;
739 }

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::Expansion3D, and Nektar::LocalRegions::Expansion2D.

Definition at line 775 of file Expansion.cpp.

780 {
781  boost::ignore_unused(dir, inarray, EdgeExp, coeffs, outarray);
782  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
783 }

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

431 {
432  const int nqtot = GetTotPoints();
433 
434  if (m_metrics.count(eMetricQuadrature) == 0)
435  {
437  }
438 
439  Vmath::Vdiv(nqtot, inarray, 1, m_metrics[eMetricQuadrature], 1, outarray,
440  1);
441 }
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.cpp:284

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

Referenced by DivideByQuadratureMetric().

◆ 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

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

Definition at line 741 of file Expansion.cpp.

745 {
746  boost::ignore_unused(data, nummodes, nmodes_offset, coeffs, fromType);
747  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
748 }

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

Referenced by ExtractDataToCoeffs().

◆ v_GenTraceExp()

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

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

Definition at line 854 of file Expansion.cpp.

855 {
856  boost::ignore_unused(traceid, exp);
858  "Method does not exist for this shape or library");
859 }

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

Referenced by GetTraceExp().

◆ v_GetCoords()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 524 of file Expansion.cpp.

527 {
528  ASSERTL1(m_geom, "m_geom not defined");
529 
530  // get physical points defined in Geom
531  m_geom->FillGeom();
532 
533  const int expDim = m_base.size();
534  int nqGeom = 1;
535  bool doCopy = true;
536 
537  Array<OneD, LibUtilities::BasisSharedPtr> CBasis(expDim);
538  Array<OneD, Array<OneD, NekDouble>> tmp(3);
539 
540  for (int i = 0; i < expDim; ++i)
541  {
542  CBasis[i] = m_geom->GetXmap()->GetBasis(i);
543  nqGeom *= CBasis[i]->GetNumPoints();
544  doCopy = doCopy &&
545  m_base[i]->GetBasisKey().SamePoints(CBasis[i]->GetBasisKey());
546  }
547 
548  tmp[0] = coords_0;
549  tmp[1] = coords_1;
550  tmp[2] = coords_2;
551 
552  if (doCopy)
553  {
554  for (int i = 0; i < m_geom->GetCoordim(); ++i)
555  {
556  m_geom->GetXmap()->BwdTrans(m_geom->GetCoeffs(i), tmp[i]);
557  }
558  }
559  else
560  {
561  for (int i = 0; i < m_geom->GetCoordim(); ++i)
562  {
563  Array<OneD, NekDouble> tmpGeom(nqGeom);
564  m_geom->GetXmap()->BwdTrans(m_geom->GetCoeffs(i), tmpGeom);
565 
566  switch (expDim)
567  {
568  case 1:
569  {
571  CBasis[0]->GetPointsKey(), &tmpGeom[0],
572  m_base[0]->GetPointsKey(), &tmp[i][0]);
573  break;
574  }
575  case 2:
576  {
578  CBasis[0]->GetPointsKey(), CBasis[1]->GetPointsKey(),
579  &tmpGeom[0], m_base[0]->GetPointsKey(),
580  m_base[1]->GetPointsKey(), &tmp[i][0]);
581  break;
582  }
583  case 3:
584  {
586  CBasis[0]->GetPointsKey(), CBasis[1]->GetPointsKey(),
587  CBasis[2]->GetPointsKey(), &tmpGeom[0],
588  m_base[0]->GetPointsKey(), m_base[1]->GetPointsKey(),
589  m_base[2]->GetPointsKey(), &tmp[i][0]);
590  break;
591  }
592  }
593  }
594  }
595 }
void Interp1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
this function interpolates a 1D function evaluated at the quadrature points of the basis fbasis0 to ...
Definition: Interp.cpp:52
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:167
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:106

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::TetExp::v_GetCoords(), Nektar::LocalRegions::TriExp::v_GetCoords(), and Nektar::LocalRegions::SegExp::v_GetCoords().

◆ v_GetLocMatrix()

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

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

Definition at line 264 of file Expansion.cpp.

266 {
267  boost::ignore_unused(mkey);
268  NEKERROR(ErrorUtil::efatal, "This function is only valid for LocalRegions");
270 }
static DNekScalMatSharedPtr NullDNekScalMatSharedPtr
Definition: NekTypeDefs.hpp:84

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

Referenced by GetLocMatrix().

◆ v_GetMF()

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

Definition at line 627 of file Expansion.cpp.

629 {
630  int coordim = GetCoordim();
631 
632  int nquad0, nquad1, nquad2;
633  int nqtot = 1;
634  switch (shapedim)
635  {
636  case 2:
637  {
638  nquad0 = m_base[0]->GetNumPoints();
639  nquad1 = m_base[1]->GetNumPoints();
640  nqtot = nquad0 * nquad1;
641  break;
642  }
643  case 3:
644  {
645  nquad0 = m_base[0]->GetNumPoints();
646  nquad1 = m_base[1]->GetNumPoints();
647  nquad2 = m_base[2]->GetNumPoints();
648  nqtot = nquad0 * nquad1 * nquad2;
649  break;
650  }
651  default:
652  break;
653  }
654 
655  StdRegions::VarCoeffType MMFCoeffs[15] = {
664 
665  Array<OneD, NekDouble> MF(coordim * nqtot);
666  Array<OneD, NekDouble> tmp(nqtot);
667  for (int k = 0; k < coordim; k++)
668  {
669  StdRegions::VarCoeffMap::const_iterator MFdir =
670  varcoeffs.find(MMFCoeffs[5 * dir + k]);
671  tmp = MFdir->second;
672 
673  Vmath::Vcopy(nqtot, &tmp[0], 1, &MF[k * nqtot], 1);
674  }
675 
676  return MF;
677 }
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255

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

◆ v_GetMFDiv()

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

◆ v_GetMFMag()

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

◆ v_GetPhysNormals()

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

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 867 of file Expansion.cpp.

868 {
869  NEKERROR(ErrorUtil::efatal, "This function is not valid for this class");
870  return NullNekDouble1DArray;
871 }

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

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

Definition at line 838 of file Expansion.cpp.

839 {
840  boost::ignore_unused(edge, outarray);
842  "Method does not exist for this shape or library");
843 }

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

Referenced by GetTracePhysMap().

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

Definition at line 828 of file Expansion.cpp.

832 {
833  boost::ignore_unused(trace, TraceExp, inarray, outarray, orient);
835  "Method does not exist for this shape or library");
836 }

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

Definition at line 820 of file Expansion.cpp.

822 {
823  boost::ignore_unused(trace, outarray);
825  "Method does not exist for this shape or library");
826 }

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 413 of file Expansion.cpp.

416 {
417  const int nqtot = GetTotPoints();
418 
419  if (m_metrics.count(eMetricQuadrature) == 0)
420  {
422  }
423 
424  Vmath::Vmul(nqtot, m_metrics[eMetricQuadrature], 1, inarray, 1, outarray,
425  1);
426 }
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:209

References 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::TriExp, Nektar::LocalRegions::QuadExp, Nektar::LocalRegions::HexExp, Nektar::LocalRegions::Expansion1D, Nektar::LocalRegions::TetExp, Nektar::LocalRegions::PyrExp, and Nektar::LocalRegions::PrismExp.

Definition at line 795 of file Expansion.cpp.

799 {
800  boost::ignore_unused(factors, d0factors, d1factors);
802  "This function is only valid for "
803  "shape expansion in LocalRegions, not parant class");
804 }

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::Expansion3D, Nektar::LocalRegions::Expansion2D, and Nektar::LocalRegions::Expansion1D.

Definition at line 845 of file Expansion.cpp.

848 {
849  boost::ignore_unused(orient, idmap, nq0, nq1);
851  "Method does not exist for this shape or library");
852 }

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 812 of file Expansion.cpp.

815 {
816  boost::ignore_unused(dir, inarray, outarray);
817  NEKERROR(ErrorUtil::efatal, "This function is not defined for this shape");
818 }

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

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

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

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

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::Expansion3D, Nektar::LocalRegions::Expansion2D, and Nektar::LocalRegions::Expansion1D.

Definition at line 901 of file Expansion.cpp.

902 {
903  boost::ignore_unused(traceid, h, p);
904  NEKERROR(ErrorUtil::efatal, "This method has not been defined");
905 }

References Nektar::ErrorUtil::efatal, NEKERROR, and CellMLToNektar.cellml_metadata::p.

Referenced by TraceNormLen().

◆ v_VectorFlux()

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

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

Definition at line 785 of file Expansion.cpp.

787 {
788  boost::ignore_unused(vec);
790  "This function is only valid for "
791  "shape expansion in LocalRegions, not parant class");
792  return 0.0;
793 }

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

145 {
146  return v_VectorFlux(vec);
147 }
virtual NekDouble v_VectorFlux(const Array< OneD, Array< OneD, NekDouble >> &vec)
Definition: Expansion.cpp:785

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

Referenced by GetRightAdjacentElementExp(), and SetAdjacentElementExp().

◆ m_elementTraceLeft

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

Definition at line 278 of file Expansion.h.

Referenced by GetLeftAdjacentElementTrace(), and SetAdjacentElementExp().

◆ m_elementTraceRight

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

Definition at line 279 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 272 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::HexExp::v_GetCoord(), Nektar::LocalRegions::PrismExp::v_GetCoord(), Nektar::LocalRegions::PyrExp::v_GetCoord(), Nektar::LocalRegions::QuadExp::v_GetCoord(), Nektar::LocalRegions::TetExp::v_GetCoord(), Nektar::LocalRegions::TriExp::v_GetCoord(), Nektar::LocalRegions::SegExp::v_GetCoord(), Nektar::LocalRegions::PrismExp::v_GetCoordim(), Nektar::LocalRegions::PyrExp::v_GetCoordim(), Nektar::LocalRegions::QuadExp::v_GetCoordim(), Nektar::LocalRegions::TetExp::v_GetCoordim(), Nektar::LocalRegions::TriExp::v_GetCoordim(), Nektar::LocalRegions::SegExp::v_GetCoordim(), Nektar::LocalRegions::PointExp::v_GetCoords(), v_GetCoords(), Nektar::LocalRegions::PrismExp::v_GetSimplexEquiSpacedConnectivity(), Nektar::LocalRegions::QuadExp::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::TriExp::v_PhysEvaluate(), Nektar::LocalRegions::SegExp::v_PhysEvaluate(), Nektar::LocalRegions::HexExp::v_PhysEvaluate(), and Nektar::LocalRegions::TetExp::v_PhysEvaluate().

◆ m_indexMapManager

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

Definition at line 269 of file Expansion.h.

Referenced by GetIndexMap().

◆ m_metricinfo

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

Definition at line 273 of file Expansion.h.

Referenced by ComputeGmatcdotMF(), ComputeQuadratureMetric(), Nektar::LocalRegions::Expansion2D::CreateMatrix(), Nektar::LocalRegions::SegExp::CreateMatrix(), CreateStaticCondMatrix(), Expansion(), GetMetricInfo(), Nektar::LocalRegions::NodalTriExp::Integral(), Nektar::LocalRegions::HexExp::IProductWRTDirectionalDerivBase_SumFac(), 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_GetMetricInfo(), 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::TetExp::v_Integral(), Nektar::LocalRegions::TriExp::v_Integral(), Nektar::LocalRegions::SegExp::v_Integral(), Nektar::LocalRegions::SegExp::v_IProductWRTBase(), Nektar::LocalRegions::SegExp::v_IProductWRTDerivBase(), Nektar::LocalRegions::TriExp::v_IProductWRTDirectionalDerivBase_SumFac(), Nektar::LocalRegions::SegExp::v_LaplacianMatrixOp(), Nektar::LocalRegions::PrismExp::v_LaplacianMatrixOp_MatFree_Kernel(), Nektar::LocalRegions::SegExp::v_MetricInfoType(), Nektar::LocalRegions::PrismExp::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