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

#include <TetExp.h>

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

Public Member Functions

 TetExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc, const SpatialDomains::TetGeomSharedPtr &geom)
 Constructor using BasisKey class for quadrature points and order definition. More...
 
 TetExp (const TetExp &T)
 Copy Constructor. More...
 
 ~TetExp ()
 Destructor. More...
 
- Public Member Functions inherited from Nektar::StdRegions::StdTetExp
 StdTetExp ()
 
 StdTetExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
 
 StdTetExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc, NekDouble *coeffs, NekDouble *phys)
 
 StdTetExp (const StdTetExp &T)
 
 ~StdTetExp ()
 
LibUtilities::ShapeType DetShapeType () const
 
NekDouble PhysEvaluate3D (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 Single Point Evaluation. More...
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion3D
 StdExpansion3D ()
 
 StdExpansion3D (int numcoeffs, const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
 
 StdExpansion3D (const StdExpansion3D &T)
 
virtual ~StdExpansion3D ()
 
void PhysTensorDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d1, Array< OneD, NekDouble > &outarray_d2, Array< OneD, NekDouble > &outarray_d3)
 Calculate the 3D derivative in the local tensor/collapsed coordinate at the physical points. More...
 
void BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
 
void IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
 
int GetNedges () const
 return the number of edges in 3D expansion More...
 
int GetEdgeNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th edge. More...
 
void GetEdgeInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards)
 
- 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 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)
 
NekDouble Linf (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( L_\infty\) error \( |\epsilon|_\infty = \max |u - u_{exact}|\) where \( u_{exact}\) is given by the array sol. More...
 
NekDouble L2 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( L_2\) error, \( | \epsilon |_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 dx \right]^{1/2} d\xi_1 \) where \( u_{exact}\) is given by the array sol. More...
 
NekDouble H1 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( H^1\) error, \( | \epsilon |^1_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 + \nabla(u - u_{exact})\cdot\nabla(u - u_{exact})\cdot dx \right]^{1/2} d\xi_1 \) where \( u_{exact}\) is given by the array sol. More...
 
const LibUtilities::PointsKeyVector GetPointsKeys () const
 
DNekMatSharedPtr BuildInverseTransformationMatrix (const DNekScalMatSharedPtr &m_transformationmatrix)
 
void PhysInterpToSimplexEquiSpaced (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int npset=-1)
 This function performs an interpolation from the physical space points provided at input into an array of equispaced points which are not the collapsed coordinate. So for a tetrahedron you will only get a tetrahedral number of values. More...
 
void GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true)
 This function provides the connectivity of local simplices (triangles or tets) to connect the equispaced data points provided by PhysInterpToSimplexEquiSpaced. More...
 
void EquiSpacedToCoeffs (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs a projection/interpolation from the equispaced points sometimes used in post-processing onto the coefficient space. More...
 
template<class T >
std::shared_ptr< T > as ()
 
void IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
void GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 
- Public Member Functions inherited from Nektar::LocalRegions::Expansion3D
 Expansion3D (SpatialDomains::Geometry3DSharedPtr pGeom)
 
virtual ~Expansion3D ()
 
void SetTraceToGeomOrientation (Array< OneD, NekDouble > &inout)
 Align trace orientation with the geometry orientation. More...
 
void SetFaceToGeomOrientation (const int face, Array< OneD, NekDouble > &inout)
 Align face orientation with the geometry orientation. More...
 
void AddHDGHelmholtzFaceTerms (const NekDouble tau, const int edge, Array< OneD, NekDouble > &facePhys, const StdRegions::VarCoeffMap &dirForcing, Array< OneD, NekDouble > &outarray)
 
void AddNormTraceInt (const int dir, Array< OneD, ExpansionSharedPtr > &FaceExp, Array< OneD, Array< OneD, NekDouble >> &faceCoeffs, Array< OneD, NekDouble > &outarray)
 
void AddNormTraceInt (const int dir, Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &FaceExp, Array< OneD, NekDouble > &outarray, const StdRegions::VarCoeffMap &varcoeffs)
 
void AddFaceBoundaryInt (const int face, ExpansionSharedPtr &FaceExp, Array< OneD, NekDouble > &facePhys, Array< OneD, NekDouble > &outarray, const StdRegions::VarCoeffMap &varcoeffs=StdRegions::NullVarCoeffMap)
 
SpatialDomains::Geometry3DSharedPtr GetGeom3D () const
 
void v_ReOrientTracePhysMap (const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1)
 
void v_NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble >> &Fvec, Array< OneD, NekDouble > &outarray)
 
Array< OneD, unsigned int > GetEdgeInverseBoundaryMap (int eid)
 
Array< OneD, unsigned int > GetTraceInverseBoundaryMap (int fid, StdRegions::Orientation faceOrient=StdRegions::eNoOrientation, int P1=-1, int P2=-1)
 
void GetInverseBoundaryMaps (Array< OneD, unsigned int > &vmap, Array< OneD, Array< OneD, unsigned int >> &emap, Array< OneD, Array< OneD, unsigned int >> &fmap)
 
DNekScalMatSharedPtr CreateMatrix (const MatrixKey &mkey)
 
- Public Member Functions inherited from Nektar::LocalRegions::Expansion
 Expansion (SpatialDomains::GeometrySharedPtr pGeom)
 
 Expansion (const Expansion &pSrc)
 
virtual ~Expansion ()
 
void SetTraceExp (const int traceid, ExpansionSharedPtr &f)
 
ExpansionSharedPtr GetTraceExp (const int traceid)
 
DNekScalMatSharedPtr GetLocMatrix (const LocalRegions::MatrixKey &mkey)
 
DNekScalMatSharedPtr GetLocMatrix (const StdRegions::MatrixType mtype, const StdRegions::ConstFactorMap &factors=StdRegions::NullConstFactorMap, const StdRegions::VarCoeffMap &varcoeffs=StdRegions::NullVarCoeffMap)
 
SpatialDomains::GeometrySharedPtr GetGeom () const
 
void Reset ()
 
IndexMapValuesSharedPtr CreateIndexMap (const IndexMapKey &ikey)
 
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)
 

Protected Member Functions

virtual NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray)
 Integrate the physical point list inarray over region. More...
 
virtual void v_PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
 Differentiate inarray in the three coordinate directions. More...
 
virtual void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Forward transform from physical quadrature space stored in inarray and evaluate the expansion coefficients and store in (this)->_coeffs. More...
 
virtual void v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Calculate the inner product of inarray with respect to the basis B=m_base0*m_base1*m_base2 and put into outarray: More...
 
virtual void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
virtual void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Calculates the inner product \( I_{pqr} = (u, \partial_{x_i} \phi_{pqr}) \). More...
 
virtual void v_AlignVectorToCollapsedDir (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
 
virtual NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
virtual NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 
virtual void v_GetCoord (const Array< OneD, const NekDouble > &Lcoords, Array< OneD, NekDouble > &coords)
 Get the coordinates "coords" at the local coordinates "Lcoords". More...
 
virtual void v_GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3)
 
virtual LibUtilities::ShapeType v_DetShapeType () const
 Return Shape of region, using ShapeType enum list. More...
 
virtual StdRegions::StdExpansionSharedPtr v_GetStdExp (void) const
 
virtual StdRegions::StdExpansionSharedPtr v_GetLinStdExp (void) const
 
virtual int v_GetCoordim ()
 
virtual void v_ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int mode_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType)
 
virtual void v_GetTracePhysMap (const int face, Array< OneD, int > &outarray)
 Returns the physical values at the quadrature points of a face. More...
 
void v_ComputeTraceNormal (const int face)
 Compute the normal of a triangular face. More...
 
virtual void v_HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
 
virtual void v_LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
 
virtual void v_LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
 
virtual void v_SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdRegions::StdMatrixKey &mkey)
 
virtual DNekMatSharedPtr v_GenMatrix (const StdRegions::StdMatrixKey &mkey)
 
virtual DNekMatSharedPtr v_CreateStdMatrix (const StdRegions::StdMatrixKey &mkey)
 
virtual DNekScalMatSharedPtr v_GetLocMatrix (const MatrixKey &mkey)
 
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const MatrixKey &mkey)
 
void v_DropLocStaticCondMatrix (const MatrixKey &mkey)
 
void SetUpInverseTransformationMatrix (const DNekMatSharedPtr &m_transformationmatrix, DNekMatSharedPtr m_inversetransformationmatrix, DNekMatSharedPtr m_inversetransposedtransformationmatrix)
 
void v_ComputeConditionNumberOfMatrix (const DNekScalMatSharedPtr &mat)
 
virtual void v_ComputeLaplacianMetric ()
 
virtual void v_NormalTraceDerivFactors (Array< OneD, Array< OneD, NekDouble >> &d0factors, Array< OneD, Array< OneD, NekDouble >> &d1factors, Array< OneD, Array< OneD, NekDouble >> &d2factors)
 : This method gets all of the factors which are required as part of the Gradient Jump Penalty stabilisation and involves the product of the normal and geometric factors along the element trace. More...
 
- Protected Member Functions inherited from Nektar::StdRegions::StdTetExp
virtual void v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
 
virtual void v_StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
 
virtual void v_IProductWRTBase_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
 
virtual void v_IProductWRTDerivBase_MatOp (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 
virtual void v_LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi)
 
virtual void v_FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 
NekDouble v_PhysEvaluateBasis (const Array< OneD, const NekDouble > &coords, int mode) final
 
virtual void v_GetTraceNumModes (const int fid, int &numModes0, int &numModes1, Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2)
 
virtual int v_GetNverts () const
 
virtual int v_GetNedges () const
 
virtual int v_GetNtraces () const
 
virtual int v_NumBndryCoeffs () const
 
virtual int v_NumDGBndryCoeffs () const
 
virtual int v_GetTraceNcoeffs (const int i) const
 
virtual int v_GetTotalTraceIntNcoeffs () const
 
virtual int v_GetTraceIntNcoeffs (const int i) const
 
virtual int v_GetTraceNumPoints (const int i) const
 
virtual int v_GetEdgeNcoeffs (const int i) const
 
virtual LibUtilities::PointsKey v_GetTracePointsKey (const int i, const int j) const
 
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
virtual const LibUtilities::BasisKey v_GetTraceBasisKey (const int i, const int k) const
 
virtual bool v_IsBoundaryInteriorExpansion ()
 
virtual int v_GetVertexMap (int localVertexId, bool useCoeffPacking=false)
 
virtual void v_GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
virtual void v_GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
virtual void v_GetTraceCoeffMap (const unsigned int fid, Array< OneD, unsigned int > &maparray)
 
virtual void v_GetElmtTraceToTraceMap (const unsigned int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
 
virtual void v_GetEdgeInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2)
 
virtual void v_GetTraceInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2)
 
virtual void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true)
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion3D
virtual NekDouble v_PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 
virtual void v_LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
 
virtual void v_HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
 
virtual void v_GetTraceToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient, int P, int Q)
 
virtual void v_GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 
- 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)
 
template<int DIR>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals)
 This function performs the barycentric interpolation of the polynomial stored in coord at a point physvals using barycentric interpolation weights in direction. More...
 
template<int DIR>
NekDouble BaryEvaluateBasis (const NekDouble &coord, const int &mode)
 
- Protected Member Functions inherited from Nektar::LocalRegions::Expansion3D
virtual void v_DGDeriv (const int dir, const Array< OneD, const NekDouble > &incoeffs, Array< OneD, ExpansionSharedPtr > &FaceExp, Array< OneD, Array< OneD, NekDouble >> &faceCoeffs, Array< OneD, NekDouble > &out_d)
 Evaluate coefficients of weak deriviative in the direction dir given the input coefficicents incoeffs and the imposed boundary values in EdgeExp (which will have its phys space updated). More...
 
virtual void v_AddFaceNormBoundaryInt (const int face, const ExpansionSharedPtr &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddRobinMassMatrix (const int face, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
virtual StdRegions::Orientation v_GetTraceOrient (int face)
 
virtual void v_GetTracePhysVals (const int face, const StdRegions::StdExpansionSharedPtr &FaceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient)
 Extract the physical values along face face from inarray into outarray following the local face orientation and point distribution defined by defined in FaceExp. More...
 
virtual void v_GenTraceExp (const int traceid, ExpansionSharedPtr &exp)
 
void GetPhysFaceVarCoeffsFromElement (const int face, ExpansionSharedPtr &FaceExp, const Array< OneD, const NekDouble > &varcoeff, Array< OneD, NekDouble > &outarray)
 
virtual Array< OneD, NekDoublev_GetnFacecdotMF (const int dir, const int face, ExpansionSharedPtr &FaceExp_f, const Array< OneD, const Array< OneD, NekDouble >> &normals, const StdRegions::VarCoeffMap &varcoeffs)
 
virtual DNekMatSharedPtr v_BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
 
virtual DNekMatSharedPtr v_BuildInverseTransformationMatrix (const DNekScalMatSharedPtr &transformationmatrix)
 Build inverse and inverse transposed transformation matrix: \(\mathbf{R^{-1}}\) and \(\mathbf{R^{-T}}\). More...
 
virtual DNekMatSharedPtr v_BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd)
 
virtual void v_TraceNormLen (const int traceid, NekDouble &h, NekDouble &p)
 
- Protected Member Functions inherited from Nektar::LocalRegions::Expansion
void ComputeLaplacianMetric ()
 
void ComputeQuadratureMetric ()
 
void ComputeGmatcdotMF (const Array< TwoD, const NekDouble > &df, const Array< OneD, const NekDouble > &direction, Array< OneD, Array< OneD, NekDouble >> &dfdir)
 
virtual void v_MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_DivideByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
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 void v_AddEdgeNormBoundaryInt (const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddEdgeNormBoundaryInt (const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddFaceNormBoundaryInt (const int face, const std::shared_ptr< Expansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
virtual NekDouble v_VectorFlux (const Array< OneD, Array< OneD, NekDouble >> &vec)
 
virtual void v_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 const Array< OneD, const NekDouble > & v_GetPhysNormals (void)
 
virtual void v_SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
virtual void v_SetUpPhysNormals (const int id)
 
virtual void v_AddRobinTraceContribution (const int traceid, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs)
 

Private Member Functions

 TetExp ()
 
void GeneralMatrixOp_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
 
virtual void v_LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
 

Private Attributes

LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLessm_matrixManager
 
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLessm_staticCondMatrixManager
 

Additional Inherited Members

- 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
 
- Protected Attributes inherited from Nektar::LocalRegions::Expansion3D
std::map< int, NormalVectorm_faceNormals
 
- Protected Attributes inherited from Nektar::LocalRegions::Expansion
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...
 

Detailed Description

Defines a Tetrahedral local expansion.

Definition at line 50 of file TetExp.h.

Constructor & Destructor Documentation

◆ TetExp() [1/3]

Nektar::LocalRegions::TetExp::TetExp ( const LibUtilities::BasisKey Ba,
const LibUtilities::BasisKey Bb,
const LibUtilities::BasisKey Bc,
const SpatialDomains::TetGeomSharedPtr geom 
)

Constructor using BasisKey class for quadrature points and order definition.

Parameters
BaBasis key for first coordinate.
BbBasis key for second coordinate.
BcBasis key for third coordinate.

Definition at line 61 of file TetExp.cpp.

66  Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
67  3, Ba, Bb, Bc),
69  Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
70  Ba, Bb, Bc),
71  StdTetExp(Ba, Bb, Bc), Expansion(geom), Expansion3D(geom),
73  std::bind(&Expansion3D::CreateMatrix, this, std::placeholders::_1),
74  std::string("TetExpMatrix")),
76  this, std::placeholders::_1),
77  std::string("TetExpStaticCondMatrix"))
78 {
79 }
Expansion3D(SpatialDomains::Geometry3DSharedPtr pGeom)
Definition: Expansion3D.h:61
DNekScalMatSharedPtr CreateMatrix(const MatrixKey &mkey)
DNekScalBlkMatSharedPtr CreateStaticCondMatrix(const MatrixKey &mkey)
Definition: Expansion.cpp:272
Expansion(SpatialDomains::GeometrySharedPtr pGeom)
Definition: Expansion.cpp:47
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: TetExp.h:196
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: TetExp.h:194
StdExpansion()
Default Constructor.
int getNumberOfCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:190

◆ TetExp() [2/3]

Nektar::LocalRegions::TetExp::TetExp ( const TetExp T)

Copy Constructor.

Definition at line 84 of file TetExp.cpp.

85  : StdRegions::StdExpansion(T), StdRegions::StdExpansion3D(T),
86  StdRegions::StdTetExp(T), Expansion(T), Expansion3D(T),
87  m_matrixManager(T.m_matrixManager),
88  m_staticCondMatrixManager(T.m_staticCondMatrixManager)
89 {
90 }

◆ ~TetExp()

Nektar::LocalRegions::TetExp::~TetExp ( )

Destructor.

Definition at line 95 of file TetExp.cpp.

96 {
97 }

◆ TetExp() [3/3]

Nektar::LocalRegions::TetExp::TetExp ( )
private

Member Function Documentation

◆ GeneralMatrixOp_MatOp()

void Nektar::LocalRegions::TetExp::GeneralMatrixOp_MatOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdRegions::StdMatrixKey mkey 
)
private

Definition at line 1068 of file TetExp.cpp.

1071 {
1072  DNekScalMatSharedPtr mat = GetLocMatrix(mkey);
1073 
1074  if (inarray.get() == outarray.get())
1075  {
1076  Array<OneD, NekDouble> tmp(m_ncoeffs);
1077  Vmath::Vcopy(m_ncoeffs, inarray.get(), 1, tmp.get(), 1);
1078 
1079  Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, mat->Scale(),
1080  (mat->GetOwnedMatrix())->GetPtr().get(), m_ncoeffs,
1081  tmp.get(), 1, 0.0, outarray.get(), 1);
1082  }
1083  else
1084  {
1085  Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, mat->Scale(),
1086  (mat->GetOwnedMatrix())->GetPtr().get(), m_ncoeffs,
1087  inarray.get(), 1, 0.0, outarray.get(), 1);
1088  }
1089 }
DNekScalMatSharedPtr GetLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:88
static void Dgemv(const char &trans, const int &m, const int &n, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = A x where A[m x n].
Definition: Blas.hpp:246
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255

References Blas::Dgemv(), Nektar::LocalRegions::Expansion::GetLocMatrix(), Nektar::StdRegions::StdExpansion::m_ncoeffs, and Vmath::Vcopy().

◆ SetUpInverseTransformationMatrix()

void Nektar::LocalRegions::TetExp::SetUpInverseTransformationMatrix ( const DNekMatSharedPtr m_transformationmatrix,
DNekMatSharedPtr  m_inversetransformationmatrix,
DNekMatSharedPtr  m_inversetransposedtransformationmatrix 
)
protected

◆ v_AlignVectorToCollapsedDir()

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

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 397 of file TetExp.cpp.

400 {
401  int i, j;
402 
403  const int nquad0 = m_base[0]->GetNumPoints();
404  const int nquad1 = m_base[1]->GetNumPoints();
405  const int nquad2 = m_base[2]->GetNumPoints();
406  const int nqtot = nquad0 * nquad1 * nquad2;
407 
408  const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
409  const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
410  const Array<OneD, const NekDouble> &z2 = m_base[2]->GetZ();
411 
412  Array<OneD, NekDouble> tmp2(nqtot);
413  Array<OneD, NekDouble> tmp3(nqtot);
414 
415  const Array<TwoD, const NekDouble> &df =
416  m_metricinfo->GetDerivFactors(GetPointsKeys());
417 
418  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
419  {
420  Vmath::Vmul(nqtot, &df[3 * dir][0], 1, inarray.get(), 1, tmp2.get(), 1);
421  Vmath::Vmul(nqtot, &df[3 * dir + 1][0], 1, inarray.get(), 1, tmp3.get(),
422  1);
423  Vmath::Vmul(nqtot, &df[3 * dir + 2][0], 1, inarray.get(), 1,
424  outarray[2].get(), 1);
425  }
426  else
427  {
428  Vmath::Smul(nqtot, df[3 * dir][0], inarray.get(), 1, tmp2.get(), 1);
429  Vmath::Smul(nqtot, df[3 * dir + 1][0], inarray.get(), 1, tmp3.get(), 1);
430  Vmath::Smul(nqtot, df[3 * dir + 2][0], inarray.get(), 1,
431  outarray[2].get(), 1);
432  }
433 
434  NekDouble g0, g1, g1a, g2, g3;
435  int k, cnt;
436 
437  for (cnt = 0, k = 0; k < nquad2; ++k)
438  {
439  for (j = 0; j < nquad1; ++j)
440  {
441  g2 = 2.0 / (1.0 - z2[k]);
442  g1 = g2 / (1.0 - z1[j]);
443  g0 = 2.0 * g1;
444  g3 = (1.0 + z1[j]) * g2 * 0.5;
445 
446  for (i = 0; i < nquad0; ++i, ++cnt)
447  {
448  g1a = g1 * (1 + z0[i]);
449 
450  outarray[0][cnt] =
451  g0 * tmp2[cnt] + g1a * (tmp3[cnt] + outarray[2][cnt]);
452 
453  outarray[1][cnt] = g2 * tmp3[cnt] + g3 * outarray[2][cnt];
454  }
455  }
456  }
457 }
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:273
const LibUtilities::PointsKeyVector GetPointsKeys() const
Array< OneD, LibUtilities::BasisSharedPtr > m_base
@ eDeformed
Geometry is curved or has non-constant factors.
double NekDouble
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:209
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:248

References Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_metricinfo, Vmath::Smul(), and Vmath::Vmul().

Referenced by v_IProductWRTDerivBase().

◆ v_ComputeConditionNumberOfMatrix()

void Nektar::LocalRegions::TetExp::v_ComputeConditionNumberOfMatrix ( const DNekScalMatSharedPtr mat)
protected

◆ v_ComputeLaplacianMetric()

void Nektar::LocalRegions::TetExp::v_ComputeLaplacianMetric ( )
protectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1169 of file TetExp.cpp.

1170 {
1171  if (m_metrics.count(eMetricQuadrature) == 0)
1172  {
1174  }
1175 
1176  int i, j;
1177  const unsigned int nqtot = GetTotPoints();
1178  const unsigned int dim = 3;
1179  const MetricType m[3][3] = {
1183 
1184  for (unsigned int i = 0; i < dim; ++i)
1185  {
1186  for (unsigned int j = i; j < dim; ++j)
1187  {
1188  m_metrics[m[i][j]] = Array<OneD, NekDouble>(nqtot);
1189  }
1190  }
1191 
1192  // Define shorthand synonyms for m_metrics storage
1193  Array<OneD, NekDouble> g0(m_metrics[m[0][0]]);
1194  Array<OneD, NekDouble> g1(m_metrics[m[1][1]]);
1195  Array<OneD, NekDouble> g2(m_metrics[m[2][2]]);
1196  Array<OneD, NekDouble> g3(m_metrics[m[0][1]]);
1197  Array<OneD, NekDouble> g4(m_metrics[m[0][2]]);
1198  Array<OneD, NekDouble> g5(m_metrics[m[1][2]]);
1199 
1200  // Allocate temporary storage
1201  Array<OneD, NekDouble> alloc(7 * nqtot, 0.0);
1202  Array<OneD, NekDouble> h0(alloc); // h0
1203  Array<OneD, NekDouble> h1(alloc + 1 * nqtot); // h1
1204  Array<OneD, NekDouble> h2(alloc + 2 * nqtot); // h2
1205  Array<OneD, NekDouble> h3(alloc + 3 * nqtot); // h3
1206  Array<OneD, NekDouble> wsp4(alloc + 4 * nqtot); // wsp4
1207  Array<OneD, NekDouble> wsp5(alloc + 5 * nqtot); // wsp5
1208  Array<OneD, NekDouble> wsp6(alloc + 6 * nqtot); // wsp6
1209  // Reuse some of the storage as workspace
1210  Array<OneD, NekDouble> wsp7(alloc); // wsp7
1211  Array<OneD, NekDouble> wsp8(alloc + 1 * nqtot); // wsp8
1212  Array<OneD, NekDouble> wsp9(alloc + 2 * nqtot); // wsp9
1213 
1214  const Array<TwoD, const NekDouble> &df =
1215  m_metricinfo->GetDerivFactors(GetPointsKeys());
1216  const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
1217  const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
1218  const Array<OneD, const NekDouble> &z2 = m_base[2]->GetZ();
1219  const unsigned int nquad0 = m_base[0]->GetNumPoints();
1220  const unsigned int nquad1 = m_base[1]->GetNumPoints();
1221  const unsigned int nquad2 = m_base[2]->GetNumPoints();
1222 
1223  for (j = 0; j < nquad2; ++j)
1224  {
1225  for (i = 0; i < nquad1; ++i)
1226  {
1227  Vmath::Fill(nquad0, 4.0 / (1.0 - z1[i]) / (1.0 - z2[j]),
1228  &h0[0] + i * nquad0 + j * nquad0 * nquad1, 1);
1229  Vmath::Fill(nquad0, 2.0 / (1.0 - z1[i]) / (1.0 - z2[j]),
1230  &h1[0] + i * nquad0 + j * nquad0 * nquad1, 1);
1231  Vmath::Fill(nquad0, 2.0 / (1.0 - z2[j]),
1232  &h2[0] + i * nquad0 + j * nquad0 * nquad1, 1);
1233  Vmath::Fill(nquad0, (1.0 + z1[i]) / (1.0 - z2[j]),
1234  &h3[0] + i * nquad0 + j * nquad0 * nquad1, 1);
1235  }
1236  }
1237  for (i = 0; i < nquad0; i++)
1238  {
1239  Blas::Dscal(nquad1 * nquad2, 1 + z0[i], &h1[0] + i, nquad0);
1240  }
1241 
1242  // Step 3. Construct combined metric terms for physical space to
1243  // collapsed coordinate system.
1244  // Order of construction optimised to minimise temporary storage
1245  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1246  {
1247  // wsp4
1248  Vmath::Vadd(nqtot, &df[1][0], 1, &df[2][0], 1, &wsp4[0], 1);
1249  Vmath::Vvtvvtp(nqtot, &df[0][0], 1, &h0[0], 1, &wsp4[0], 1, &h1[0], 1,
1250  &wsp4[0], 1);
1251  // wsp5
1252  Vmath::Vadd(nqtot, &df[4][0], 1, &df[5][0], 1, &wsp5[0], 1);
1253  Vmath::Vvtvvtp(nqtot, &df[3][0], 1, &h0[0], 1, &wsp5[0], 1, &h1[0], 1,
1254  &wsp5[0], 1);
1255  // wsp6
1256  Vmath::Vadd(nqtot, &df[7][0], 1, &df[8][0], 1, &wsp6[0], 1);
1257  Vmath::Vvtvvtp(nqtot, &df[6][0], 1, &h0[0], 1, &wsp6[0], 1, &h1[0], 1,
1258  &wsp6[0], 1);
1259 
1260  // g0
1261  Vmath::Vvtvvtp(nqtot, &wsp4[0], 1, &wsp4[0], 1, &wsp5[0], 1, &wsp5[0],
1262  1, &g0[0], 1);
1263  Vmath::Vvtvp(nqtot, &wsp6[0], 1, &wsp6[0], 1, &g0[0], 1, &g0[0], 1);
1264 
1265  // g4
1266  Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &wsp4[0], 1, &df[5][0], 1, &wsp5[0],
1267  1, &g4[0], 1);
1268  Vmath::Vvtvp(nqtot, &df[8][0], 1, &wsp6[0], 1, &g4[0], 1, &g4[0], 1);
1269 
1270  // overwrite h0, h1, h2
1271  // wsp7 (h2f1 + h3f2)
1272  Vmath::Vvtvvtp(nqtot, &df[1][0], 1, &h2[0], 1, &df[2][0], 1, &h3[0], 1,
1273  &wsp7[0], 1);
1274  // wsp8 (h2f4 + h3f5)
1275  Vmath::Vvtvvtp(nqtot, &df[4][0], 1, &h2[0], 1, &df[5][0], 1, &h3[0], 1,
1276  &wsp8[0], 1);
1277  // wsp9 (h2f7 + h3f8)
1278  Vmath::Vvtvvtp(nqtot, &df[7][0], 1, &h2[0], 1, &df[8][0], 1, &h3[0], 1,
1279  &wsp9[0], 1);
1280 
1281  // g3
1282  Vmath::Vvtvvtp(nqtot, &wsp4[0], 1, &wsp7[0], 1, &wsp5[0], 1, &wsp8[0],
1283  1, &g3[0], 1);
1284  Vmath::Vvtvp(nqtot, &wsp6[0], 1, &wsp9[0], 1, &g3[0], 1, &g3[0], 1);
1285 
1286  // overwrite wsp4, wsp5, wsp6
1287  // g1
1288  Vmath::Vvtvvtp(nqtot, &wsp7[0], 1, &wsp7[0], 1, &wsp8[0], 1, &wsp8[0],
1289  1, &g1[0], 1);
1290  Vmath::Vvtvp(nqtot, &wsp9[0], 1, &wsp9[0], 1, &g1[0], 1, &g1[0], 1);
1291 
1292  // g5
1293  Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &wsp7[0], 1, &df[5][0], 1, &wsp8[0],
1294  1, &g5[0], 1);
1295  Vmath::Vvtvp(nqtot, &df[8][0], 1, &wsp9[0], 1, &g5[0], 1, &g5[0], 1);
1296 
1297  // g2
1298  Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &df[2][0], 1, &df[5][0], 1,
1299  &df[5][0], 1, &g2[0], 1);
1300  Vmath::Vvtvp(nqtot, &df[8][0], 1, &df[8][0], 1, &g2[0], 1, &g2[0], 1);
1301  }
1302  else
1303  {
1304  // wsp4
1305  Vmath::Svtsvtp(nqtot, df[0][0], &h0[0], 1, df[1][0] + df[2][0], &h1[0],
1306  1, &wsp4[0], 1);
1307  // wsp5
1308  Vmath::Svtsvtp(nqtot, df[3][0], &h0[0], 1, df[4][0] + df[5][0], &h1[0],
1309  1, &wsp5[0], 1);
1310  // wsp6
1311  Vmath::Svtsvtp(nqtot, df[6][0], &h0[0], 1, df[7][0] + df[8][0], &h1[0],
1312  1, &wsp6[0], 1);
1313 
1314  // g0
1315  Vmath::Vvtvvtp(nqtot, &wsp4[0], 1, &wsp4[0], 1, &wsp5[0], 1, &wsp5[0],
1316  1, &g0[0], 1);
1317  Vmath::Vvtvp(nqtot, &wsp6[0], 1, &wsp6[0], 1, &g0[0], 1, &g0[0], 1);
1318 
1319  // g4
1320  Vmath::Svtsvtp(nqtot, df[2][0], &wsp4[0], 1, df[5][0], &wsp5[0], 1,
1321  &g4[0], 1);
1322  Vmath::Svtvp(nqtot, df[8][0], &wsp6[0], 1, &g4[0], 1, &g4[0], 1);
1323 
1324  // overwrite h0, h1, h2
1325  // wsp7 (h2f1 + h3f2)
1326  Vmath::Svtsvtp(nqtot, df[1][0], &h2[0], 1, df[2][0], &h3[0], 1,
1327  &wsp7[0], 1);
1328  // wsp8 (h2f4 + h3f5)
1329  Vmath::Svtsvtp(nqtot, df[4][0], &h2[0], 1, df[5][0], &h3[0], 1,
1330  &wsp8[0], 1);
1331  // wsp9 (h2f7 + h3f8)
1332  Vmath::Svtsvtp(nqtot, df[7][0], &h2[0], 1, df[8][0], &h3[0], 1,
1333  &wsp9[0], 1);
1334 
1335  // g3
1336  Vmath::Vvtvvtp(nqtot, &wsp4[0], 1, &wsp7[0], 1, &wsp5[0], 1, &wsp8[0],
1337  1, &g3[0], 1);
1338  Vmath::Vvtvp(nqtot, &wsp6[0], 1, &wsp9[0], 1, &g3[0], 1, &g3[0], 1);
1339 
1340  // overwrite wsp4, wsp5, wsp6
1341  // g1
1342  Vmath::Vvtvvtp(nqtot, &wsp7[0], 1, &wsp7[0], 1, &wsp8[0], 1, &wsp8[0],
1343  1, &g1[0], 1);
1344  Vmath::Vvtvp(nqtot, &wsp9[0], 1, &wsp9[0], 1, &g1[0], 1, &g1[0], 1);
1345 
1346  // g5
1347  Vmath::Svtsvtp(nqtot, df[2][0], &wsp7[0], 1, df[5][0], &wsp8[0], 1,
1348  &g5[0], 1);
1349  Vmath::Svtvp(nqtot, df[8][0], &wsp9[0], 1, &g5[0], 1, &g5[0], 1);
1350 
1351  // g2
1352  Vmath::Fill(nqtot,
1353  df[2][0] * df[2][0] + df[5][0] * df[5][0] +
1354  df[8][0] * df[8][0],
1355  &g2[0], 1);
1356  }
1357 
1358  for (unsigned int i = 0; i < dim; ++i)
1359  {
1360  for (unsigned int j = i; j < dim; ++j)
1361  {
1362  MultiplyByQuadratureMetric(m_metrics[m[i][j]], m_metrics[m[i][j]]);
1363  }
1364  }
1365 }
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:731
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition: Blas.hpp:168
void Svtsvtp(int n, const T alpha, const T *x, int incx, const T beta, const T *y, int incy, T *z, int incz)
vvtvvtp (scalar times vector plus scalar times vector):
Definition: Vmath.cpp:751
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
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:359
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45
void Vvtvvtp(int n, const T *v, int incv, const T *w, int incw, const T *x, int incx, const T *y, int incy, T *z, int incz)
vvtvvtp (vector times vector plus vector times vector):
Definition: Vmath.cpp:692

◆ v_ComputeTraceNormal()

void Nektar::LocalRegions::TetExp::v_ComputeTraceNormal ( const int  face)
protectedvirtual

Compute the normal of a triangular face.

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 708 of file TetExp.cpp.

709 {
710  int i;
711  const SpatialDomains::GeomFactorsSharedPtr &geomFactors =
712  GetGeom()->GetMetricInfo();
713 
715  for (int i = 0; i < ptsKeys.size(); ++i)
716  {
717  // Need at least 2 points for computing normals
718  if (ptsKeys[i].GetNumPoints() == 1)
719  {
720  LibUtilities::PointsKey pKey(2, ptsKeys[i].GetPointsType());
721  ptsKeys[i] = pKey;
722  }
723  }
724 
725  SpatialDomains::GeomType type = geomFactors->GetGtype();
726  const Array<TwoD, const NekDouble> &df =
727  geomFactors->GetDerivFactors(ptsKeys);
728  const Array<OneD, const NekDouble> &jac = geomFactors->GetJac(ptsKeys);
729 
730  LibUtilities::BasisKey tobasis0 = GetTraceBasisKey(face, 0);
731  LibUtilities::BasisKey tobasis1 = GetTraceBasisKey(face, 1);
732 
733  // number of face quadrature points
734  int nq_face = tobasis0.GetNumPoints() * tobasis1.GetNumPoints();
735 
736  int vCoordDim = GetCoordim();
737 
738  m_traceNormals[face] = Array<OneD, Array<OneD, NekDouble>>(vCoordDim);
739  Array<OneD, Array<OneD, NekDouble>> &normal = m_traceNormals[face];
740  for (i = 0; i < vCoordDim; ++i)
741  {
742  normal[i] = Array<OneD, NekDouble>(nq_face);
743  }
744 
745  size_t nqb = nq_face;
746  size_t nbnd = face;
747  m_elmtBndNormDirElmtLen[nbnd] = Array<OneD, NekDouble>{nqb, 0.0};
748  Array<OneD, NekDouble> &length = m_elmtBndNormDirElmtLen[nbnd];
749 
750  // Regular geometry case
751  if (type == SpatialDomains::eRegular ||
753  {
754  NekDouble fac;
755 
756  // Set up normals
757  switch (face)
758  {
759  case 0:
760  {
761  for (i = 0; i < vCoordDim; ++i)
762  {
763  normal[i][0] = -df[3 * i + 2][0];
764  }
765 
766  break;
767  }
768  case 1:
769  {
770  for (i = 0; i < vCoordDim; ++i)
771  {
772  normal[i][0] = -df[3 * i + 1][0];
773  }
774 
775  break;
776  }
777  case 2:
778  {
779  for (i = 0; i < vCoordDim; ++i)
780  {
781  normal[i][0] =
782  df[3 * i][0] + df[3 * i + 1][0] + df[3 * i + 2][0];
783  }
784 
785  break;
786  }
787  case 3:
788  {
789  for (i = 0; i < vCoordDim; ++i)
790  {
791  normal[i][0] = -df[3 * i][0];
792  }
793  break;
794  }
795  default:
796  ASSERTL0(false, "face is out of range (edge < 3)");
797  }
798 
799  // normalise
800  fac = 0.0;
801  for (i = 0; i < vCoordDim; ++i)
802  {
803  fac += normal[i][0] * normal[i][0];
804  }
805  fac = 1.0 / sqrt(fac);
806  Vmath::Fill(nqb, fac, length, 1);
807 
808  for (i = 0; i < vCoordDim; ++i)
809  {
810  Vmath::Fill(nq_face, fac * normal[i][0], normal[i], 1);
811  }
812  }
813  else
814  {
815  // Set up deformed normals
816  int j, k;
817 
818  int nq0 = ptsKeys[0].GetNumPoints();
819  int nq1 = ptsKeys[1].GetNumPoints();
820  int nq2 = ptsKeys[2].GetNumPoints();
821  int nqtot;
822  int nq01 = nq0 * nq1;
823 
824  // number of elemental quad points
825  if (face == 0)
826  {
827  nqtot = nq01;
828  }
829  else if (face == 1)
830  {
831  nqtot = nq0 * nq2;
832  }
833  else
834  {
835  nqtot = nq1 * nq2;
836  }
837 
838  LibUtilities::PointsKey points0;
839  LibUtilities::PointsKey points1;
840 
841  Array<OneD, NekDouble> faceJac(nqtot);
842  Array<OneD, NekDouble> normals(vCoordDim * nqtot, 0.0);
843 
844  // Extract Jacobian along face and recover local derivates
845  // (dx/dr) for polynomial interpolation by multiplying m_gmat by
846  // jacobian
847  switch (face)
848  {
849  case 0:
850  {
851  for (j = 0; j < nq01; ++j)
852  {
853  normals[j] = -df[2][j] * jac[j];
854  normals[nqtot + j] = -df[5][j] * jac[j];
855  normals[2 * nqtot + j] = -df[8][j] * jac[j];
856  faceJac[j] = jac[j];
857  }
858 
859  points0 = ptsKeys[0];
860  points1 = ptsKeys[1];
861  break;
862  }
863 
864  case 1:
865  {
866  for (j = 0; j < nq0; ++j)
867  {
868  for (k = 0; k < nq2; ++k)
869  {
870  int tmp = j + nq01 * k;
871  normals[j + k * nq0] = -df[1][tmp] * jac[tmp];
872  normals[nqtot + j + k * nq0] = -df[4][tmp] * jac[tmp];
873  normals[2 * nqtot + j + k * nq0] =
874  -df[7][tmp] * jac[tmp];
875  faceJac[j + k * nq0] = jac[tmp];
876  }
877  }
878 
879  points0 = ptsKeys[0];
880  points1 = ptsKeys[2];
881  break;
882  }
883 
884  case 2:
885  {
886  for (j = 0; j < nq1; ++j)
887  {
888  for (k = 0; k < nq2; ++k)
889  {
890  int tmp = nq0 - 1 + nq0 * j + nq01 * k;
891  normals[j + k * nq1] =
892  (df[0][tmp] + df[1][tmp] + df[2][tmp]) * jac[tmp];
893  normals[nqtot + j + k * nq1] =
894  (df[3][tmp] + df[4][tmp] + df[5][tmp]) * jac[tmp];
895  normals[2 * nqtot + j + k * nq1] =
896  (df[6][tmp] + df[7][tmp] + df[8][tmp]) * jac[tmp];
897  faceJac[j + k * nq1] = jac[tmp];
898  }
899  }
900 
901  points0 = ptsKeys[1];
902  points1 = ptsKeys[2];
903  break;
904  }
905 
906  case 3:
907  {
908  for (j = 0; j < nq1; ++j)
909  {
910  for (k = 0; k < nq2; ++k)
911  {
912  int tmp = j * nq0 + nq01 * k;
913  normals[j + k * nq1] = -df[0][tmp] * jac[tmp];
914  normals[nqtot + j + k * nq1] = -df[3][tmp] * jac[tmp];
915  normals[2 * nqtot + j + k * nq1] =
916  -df[6][tmp] * jac[tmp];
917  faceJac[j + k * nq1] = jac[tmp];
918  }
919  }
920 
921  points0 = ptsKeys[1];
922  points1 = ptsKeys[2];
923  break;
924  }
925 
926  default:
927  ASSERTL0(false, "face is out of range (face < 3)");
928  }
929 
930  Array<OneD, NekDouble> work(nq_face, 0.0);
931  // Interpolate Jacobian and invert
932  LibUtilities::Interp2D(points0, points1, faceJac,
933  tobasis0.GetPointsKey(), tobasis1.GetPointsKey(),
934  work);
935  Vmath::Sdiv(nq_face, 1.0, &work[0], 1, &work[0], 1);
936 
937  // Interpolate normal and multiply by inverse Jacobian.
938  for (i = 0; i < vCoordDim; ++i)
939  {
940  LibUtilities::Interp2D(points0, points1, &normals[i * nqtot],
941  tobasis0.GetPointsKey(),
942  tobasis1.GetPointsKey(), &normal[i][0]);
943  Vmath::Vmul(nq_face, work, 1, normal[i], 1, normal[i], 1);
944  }
945 
946  // Normalise to obtain unit normals.
947  Vmath::Zero(nq_face, work, 1);
948  for (i = 0; i < GetCoordim(); ++i)
949  {
950  Vmath::Vvtvp(nq_face, normal[i], 1, normal[i], 1, work, 1, work, 1);
951  }
952 
953  Vmath::Vsqrt(nq_face, work, 1, work, 1);
954  Vmath::Sdiv(nq_face, 1.0, work, 1, work, 1);
955 
956  Vmath::Vcopy(nqb, work, 1, length, 1);
957 
958  for (i = 0; i < GetCoordim(); ++i)
959  {
960  Vmath::Vmul(nq_face, normal[i], 1, work, 1, normal[i], 1);
961  }
962  }
963 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
std::map< int, NormalVector > m_traceNormals
Definition: Expansion.h:275
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
SpatialDomains::GeometrySharedPtr GetGeom() const
Definition: Expansion.cpp:166
const LibUtilities::BasisKey GetTraceBasisKey(const int i, int k=-1) const
This function returns the basis key belonging to the i-th trace.
Definition: StdExpansion.h:307
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:213
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:226
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
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:250
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition: GeomFactors.h:62
GeomType
Indicates the type of element geometry.
@ eRegular
Geometry is straight-sided with constant geometric factors.
@ eMovingRegular
Currently unused.
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
Definition: Vmath.cpp:534
void Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/y.
Definition: Vmath.cpp:324
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:291

References ASSERTL0, Nektar::SpatialDomains::eMovingRegular, Nektar::SpatialDomains::eRegular, Vmath::Fill(), Nektar::StdRegions::StdExpansion::GetCoordim(), Nektar::LocalRegions::Expansion::GetGeom(), Nektar::LibUtilities::BasisKey::GetNumPoints(), Nektar::StdRegions::StdExpansion::GetNumPoints(), Nektar::LibUtilities::BasisKey::GetPointsKey(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::StdRegions::StdExpansion::GetTraceBasisKey(), Nektar::LibUtilities::Interp2D(), Nektar::LocalRegions::Expansion::m_elmtBndNormDirElmtLen, Nektar::LocalRegions::Expansion::m_traceNormals, Vmath::Sdiv(), tinysimd::sqrt(), Vmath::Vcopy(), Vmath::Vmul(), Vmath::Vsqrt(), Vmath::Vvtvp(), and Vmath::Zero().

◆ v_CreateStdMatrix()

DNekMatSharedPtr Nektar::LocalRegions::TetExp::v_CreateStdMatrix ( const StdRegions::StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdTetExp.

Definition at line 1042 of file TetExp.cpp.

1043 {
1044  LibUtilities::BasisKey bkey0 = m_base[0]->GetBasisKey();
1045  LibUtilities::BasisKey bkey1 = m_base[1]->GetBasisKey();
1046  LibUtilities::BasisKey bkey2 = m_base[2]->GetBasisKey();
1048  MemoryManager<StdTetExp>::AllocateSharedPtr(bkey0, bkey1, bkey2);
1049 
1050  return tmp->GetStdMatrix(mkey);
1051 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< StdTetExp > StdTetExpSharedPtr
Definition: StdTetExp.h:247

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and Nektar::StdRegions::StdExpansion::m_base.

◆ v_DetShapeType()

LibUtilities::ShapeType Nektar::LocalRegions::TetExp::v_DetShapeType ( ) const
protectedvirtual

Return Shape of region, using ShapeType enum list.

Reimplemented from Nektar::StdRegions::StdTetExp.

Definition at line 528 of file TetExp.cpp.

529 {
531 }

References Nektar::LibUtilities::eTetrahedron.

◆ v_DropLocStaticCondMatrix()

void Nektar::LocalRegions::TetExp::v_DropLocStaticCondMatrix ( const MatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1063 of file TetExp.cpp.

1064 {
1065  m_staticCondMatrixManager.DeleteObject(mkey);
1066 }

References m_staticCondMatrixManager.

◆ v_ExtractDataToCoeffs()

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

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 558 of file TetExp.cpp.

562 {
563  boost::ignore_unused(fromType);
564 
565  int data_order0 = nummodes[mode_offset];
566  int fillorder0 = min(m_base[0]->GetNumModes(), data_order0);
567  int data_order1 = nummodes[mode_offset + 1];
568  int order1 = m_base[1]->GetNumModes();
569  int fillorder1 = min(order1, data_order1);
570  int data_order2 = nummodes[mode_offset + 2];
571  int order2 = m_base[2]->GetNumModes();
572  int fillorder2 = min(order2, data_order2);
573 
574  switch (m_base[0]->GetBasisType())
575  {
577  {
578  int i, j;
579  int cnt = 0;
580  int cnt1 = 0;
581 
583  "Extraction routine not set up for this basis");
585  "Extraction routine not set up for this basis");
586 
587  Vmath::Zero(m_ncoeffs, coeffs, 1);
588  for (j = 0; j < fillorder0; ++j)
589  {
590  for (i = 0; i < fillorder1 - j; ++i)
591  {
592  Vmath::Vcopy(fillorder2 - j - i, &data[cnt], 1,
593  &coeffs[cnt1], 1);
594  cnt += data_order2 - j - i;
595  cnt1 += order2 - j - i;
596  }
597 
598  // count out data for j iteration
599  for (i = fillorder1 - j; i < data_order1 - j; ++i)
600  {
601  cnt += data_order2 - j - i;
602  }
603 
604  for (i = fillorder1 - j; i < order1 - j; ++i)
605  {
606  cnt1 += order2 - j - i;
607  }
608  }
609  }
610  break;
611  default:
612  ASSERTL0(false, "basis is either not set up or not "
613  "hierarchicial");
614  }
615 }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:163
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:51
@ eModified_C
Principle Modified Functions .
Definition: BasisType.h:52
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:50

References ASSERTL0, ASSERTL1, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::LibUtilities::eModified_C, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Vmath::Vcopy(), and Vmath::Zero().

◆ v_FwdTrans()

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

Forward transform from physical quadrature space stored in inarray and evaluate the expansion coefficients and store in (this)->_coeffs.

Parameters
inarrayArray of physical quadrature points to be transformed.
outarrayArray of coefficients to update.

Reimplemented from Nektar::StdRegions::StdTetExp.

Definition at line 234 of file TetExp.cpp.

236 {
237  if ((m_base[0]->Collocation()) && (m_base[1]->Collocation()) &&
238  (m_base[2]->Collocation()))
239  {
240  Vmath::Vcopy(GetNcoeffs(), &inarray[0], 1, &outarray[0], 1);
241  }
242  else
243  {
244  IProductWRTBase(inarray, outarray);
245 
246  // get Mass matrix inverse
247  MatrixKey masskey(StdRegions::eInvMass, DetShapeType(), *this);
248  DNekScalMatSharedPtr matsys = m_matrixManager[masskey];
249 
250  // copy inarray in case inarray == outarray
251  DNekVec in(m_ncoeffs, outarray);
252  DNekVec out(m_ncoeffs, outarray, eWrapper);
253 
254  out = (*matsys) * in;
255  }
256 }
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:130
void IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
this function calculates the inner product of a given function f with the different modes of the expa...
Definition: StdExpansion.h:536
LibUtilities::ShapeType DetShapeType() const
Definition: StdTetExp.h:64
NekVector< NekDouble > DNekVec
Definition: NekTypeDefs.hpp:48

References Nektar::StdRegions::StdTetExp::DetShapeType(), Nektar::StdRegions::eInvMass, Nektar::eWrapper, Nektar::StdRegions::StdExpansion::GetNcoeffs(), Nektar::StdRegions::StdExpansion::IProductWRTBase(), Nektar::StdRegions::StdExpansion::m_base, m_matrixManager, Nektar::StdRegions::StdExpansion::m_ncoeffs, and Vmath::Vcopy().

◆ v_GenMatrix()

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

Reimplemented from Nektar::StdRegions::StdTetExp.

Definition at line 1020 of file TetExp.cpp.

1021 {
1022  DNekMatSharedPtr returnval;
1023 
1024  switch (mkey.GetMatrixType())
1025  {
1033  returnval = Expansion3D::v_GenMatrix(mkey);
1034  break;
1035  default:
1036  returnval = StdTetExp::v_GenMatrix(mkey);
1037  }
1038 
1039  return returnval;
1040 }
virtual DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey)
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75

References Nektar::StdRegions::eHybridDGHelmBndLam, Nektar::StdRegions::eHybridDGHelmholtz, Nektar::StdRegions::eHybridDGLamToQ0, Nektar::StdRegions::eHybridDGLamToQ1, Nektar::StdRegions::eHybridDGLamToQ2, Nektar::StdRegions::eHybridDGLamToU, Nektar::StdRegions::eInvLaplacianWithUnityMean, Nektar::StdRegions::StdMatrixKey::GetMatrixType(), and Nektar::LocalRegions::Expansion3D::v_GenMatrix().

◆ v_GetCoord()

void Nektar::LocalRegions::TetExp::v_GetCoord ( const Array< OneD, const NekDouble > &  Lcoords,
Array< OneD, NekDouble > &  coords 
)
protectedvirtual

Get the coordinates "coords" at the local coordinates "Lcoords".

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 497 of file TetExp.cpp.

499 {
500  int i;
501 
502  ASSERTL1(Lcoords[0] <= -1.0 && Lcoords[0] >= 1.0 && Lcoords[1] <= -1.0 &&
503  Lcoords[1] >= 1.0 && Lcoords[2] <= -1.0 && Lcoords[2] >= 1.0,
504  "Local coordinates are not in region [-1,1]");
505 
506  // m_geom->FillGeom(); // TODO: implement FillGeom()
507 
508  for (i = 0; i < m_geom->GetCoordim(); ++i)
509  {
510  coords[i] = m_geom->GetCoord(i, Lcoords);
511  }
512 }
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:272

References ASSERTL1, and Nektar::LocalRegions::Expansion::m_geom.

◆ v_GetCoordim()

int Nektar::LocalRegions::TetExp::v_GetCoordim ( void  )
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 553 of file TetExp.cpp.

554 {
555  return m_geom->GetCoordim();
556 }

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

◆ v_GetCoords()

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

Reimplemented from Nektar::StdRegions::StdTetExp.

Definition at line 514 of file TetExp.cpp.

517 {
518  Expansion::v_GetCoords(coords_0, coords_1, coords_2);
519 }
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3)
Definition: Expansion.cpp:524

References Nektar::LocalRegions::Expansion::v_GetCoords().

◆ v_GetLinStdExp()

StdRegions::StdExpansionSharedPtr Nektar::LocalRegions::TetExp::v_GetLinStdExp ( void  ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 540 of file TetExp.cpp.

541 {
542  LibUtilities::BasisKey bkey0(m_base[0]->GetBasisType(), 2,
543  m_base[0]->GetPointsKey());
544  LibUtilities::BasisKey bkey1(m_base[1]->GetBasisType(), 2,
545  m_base[1]->GetPointsKey());
546  LibUtilities::BasisKey bkey2(m_base[2]->GetBasisType(), 2,
547  m_base[2]->GetPointsKey());
548 
550  bkey2);
551 }

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::StdRegions::StdExpansion::GetBasisType(), and Nektar::StdRegions::StdExpansion::m_base.

◆ v_GetLocMatrix()

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

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1053 of file TetExp.cpp.

1054 {
1055  return m_matrixManager[mkey];
1056 }

References m_matrixManager.

◆ v_GetLocStaticCondMatrix()

DNekScalBlkMatSharedPtr Nektar::LocalRegions::TetExp::v_GetLocStaticCondMatrix ( const MatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1058 of file TetExp.cpp.

1059 {
1060  return m_staticCondMatrixManager[mkey];
1061 }

References m_staticCondMatrixManager.

◆ v_GetStdExp()

StdRegions::StdExpansionSharedPtr Nektar::LocalRegions::TetExp::v_GetStdExp ( void  ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 533 of file TetExp.cpp.

534 {
536  m_base[0]->GetBasisKey(), m_base[1]->GetBasisKey(),
537  m_base[2]->GetBasisKey());
538 }

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and Nektar::StdRegions::StdExpansion::m_base.

◆ v_GetTracePhysMap()

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

Returns the physical values at the quadrature points of a face.

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 620 of file TetExp.cpp.

621 {
622  int nquad0 = m_base[0]->GetNumPoints();
623  int nquad1 = m_base[1]->GetNumPoints();
624  int nquad2 = m_base[2]->GetNumPoints();
625 
626  int nq0 = 0;
627  int nq1 = 0;
628 
629  // get forward aligned faces.
630  switch (face)
631  {
632  case 0:
633  {
634  nq0 = nquad0;
635  nq1 = nquad1;
636  if (outarray.size() != nq0 * nq1)
637  {
638  outarray = Array<OneD, int>(nq0 * nq1);
639  }
640 
641  for (int i = 0; i < nquad0 * nquad1; ++i)
642  {
643  outarray[i] = i;
644  }
645 
646  break;
647  }
648  case 1:
649  {
650  nq0 = nquad0;
651  nq1 = nquad2;
652  if (outarray.size() != nq0 * nq1)
653  {
654  outarray = Array<OneD, int>(nq0 * nq1);
655  }
656 
657  // Direction A and B positive
658  for (int k = 0; k < nquad2; k++)
659  {
660  for (int i = 0; i < nquad0; ++i)
661  {
662  outarray[k * nquad0 + i] = (nquad0 * nquad1 * k) + i;
663  }
664  }
665  break;
666  }
667  case 2:
668  {
669  nq0 = nquad1;
670  nq1 = nquad2;
671  if (outarray.size() != nq0 * nq1)
672  {
673  outarray = Array<OneD, int>(nq0 * nq1);
674  }
675 
676  // Directions A and B positive
677  for (int j = 0; j < nquad1 * nquad2; ++j)
678  {
679  outarray[j] = nquad0 - 1 + j * nquad0;
680  }
681  break;
682  }
683  case 3:
684  {
685  nq0 = nquad1;
686  nq1 = nquad2;
687  if (outarray.size() != nq0 * nq1)
688  {
689  outarray = Array<OneD, int>(nq0 * nq1);
690  }
691 
692  // Directions A and B positive
693  for (int j = 0; j < nquad1 * nquad2; ++j)
694  {
695  outarray[j] = j * nquad0;
696  }
697  }
698  break;
699  default:
700  ASSERTL0(false, "face value (> 3) is out of range");
701  break;
702  }
703 }

References ASSERTL0, and Nektar::StdRegions::StdExpansion::m_base.

◆ v_HelmholtzMatrixOp()

void Nektar::LocalRegions::TetExp::v_HelmholtzMatrixOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdRegions::StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 968 of file TetExp.cpp.

971 {
972  TetExp::v_HelmholtzMatrixOp_MatFree(inarray, outarray, mkey);
973 }
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)

References Nektar::StdRegions::StdExpansion3D::v_HelmholtzMatrixOp_MatFree().

◆ v_Integral()

NekDouble Nektar::LocalRegions::TetExp::v_Integral ( const Array< OneD, const NekDouble > &  inarray)
protectedvirtual

Integrate the physical point list inarray over region.

Parameters
inarrayDefinition of function to be returned at quadrature point of expansion.
Returns
\(\int^1_{-1}\int^1_{-1} \int^1_{-1} u(\eta_1, \eta_2, \eta_3) J[i,j,k] d \eta_1 d \eta_2 d \eta_3 \) where \(inarray[i,j,k] = u(\eta_{1i},\eta_{2j},\eta_{3k}) \) and \( J[i,j,k] \) is the Jacobian evaluated at the quadrature point.

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 113 of file TetExp.cpp.

114 {
115  int nquad0 = m_base[0]->GetNumPoints();
116  int nquad1 = m_base[1]->GetNumPoints();
117  int nquad2 = m_base[2]->GetNumPoints();
118  Array<OneD, const NekDouble> jac = m_metricinfo->GetJac(GetPointsKeys());
119  NekDouble retrunVal;
120  Array<OneD, NekDouble> tmp(nquad0 * nquad1 * nquad2);
121 
122  // multiply inarray with Jacobian
123  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
124  {
125  Vmath::Vmul(nquad0 * nquad1 * nquad2, &jac[0], 1,
126  (NekDouble *)&inarray[0], 1, &tmp[0], 1);
127  }
128  else
129  {
130  Vmath::Smul(nquad0 * nquad1 * nquad2, (NekDouble)jac[0],
131  (NekDouble *)&inarray[0], 1, &tmp[0], 1);
132  }
133 
134  // call StdTetExp version;
135  retrunVal = StdTetExp::v_Integral(tmp);
136 
137  return retrunVal;
138 }

References Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_metricinfo, Vmath::Smul(), and Vmath::Vmul().

◆ v_IProductWRTBase()

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

Calculate the inner product of inarray with respect to the basis B=m_base0*m_base1*m_base2 and put into outarray:

\( \begin{array}{rcl} I_{pqr} = (\phi_{pqr}, u)_{\delta} & = & \sum_{i=0}^{nq_0} \sum_{j=0}^{nq_1} \sum_{k=0}^{nq_2} \psi_{p}^{a} (\eta_{1i}) \psi_{pq}^{b} (\eta_{2j}) \psi_{pqr}^{c} (\eta_{3k}) w_i w_j w_k u(\eta_{1,i} \eta_{2,j} \eta_{3,k}) J_{i,j,k}\\ & = & \sum_{i=0}^{nq_0} \psi_p^a(\eta_{1,i}) \sum_{j=0}^{nq_1} \psi_{pq}^b(\eta_{2,j}) \sum_{k=0}^{nq_2} \psi_{pqr}^c u(\eta_{1i},\eta_{2j},\eta_{3k}) J_{i,j,k} \end{array} \)
where \( \phi_{pqr} (\xi_1 , \xi_2 , \xi_3) = \psi_p^a (\eta_1) \psi_{pq}^b (\eta_2) \psi_{pqr}^c (\eta_3) \) which can be implemented as
\(f_{pqr} (\xi_{3k}) = \sum_{k=0}^{nq_3} \psi_{pqr}^c u(\eta_{1i},\eta_{2j},\eta_{3k}) J_{i,j,k} = {\bf B_3 U} \)
\( g_{pq} (\xi_{3k}) = \sum_{j=0}^{nq_1} \psi_{pq}^b (\xi_{2j}) f_{pqr} (\xi_{3k}) = {\bf B_2 F} \)
\( (\phi_{pqr}, u)_{\delta} = \sum_{k=0}^{nq_0} \psi_{p}^a (\xi_{3k}) g_{pq} (\xi_{3k}) = {\bf B_1 G} \)

Reimplemented from Nektar::StdRegions::StdTetExp.

Definition at line 287 of file TetExp.cpp.

289 {
290  v_IProductWRTBase_SumFac(inarray, outarray);
291 }
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
Definition: TetExp.cpp:293

References v_IProductWRTBase_SumFac().

◆ v_IProductWRTBase_SumFac()

void Nektar::LocalRegions::TetExp::v_IProductWRTBase_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
bool  multiplybyweights = true 
)
protectedvirtual
Parameters
inarrayFunction evaluated at physical collocation points.
outarrayInner product with respect to each basis function over the element.

Reimplemented from Nektar::StdRegions::StdTetExp.

Definition at line 293 of file TetExp.cpp.

296 {
297  const int nquad0 = m_base[0]->GetNumPoints();
298  const int nquad1 = m_base[1]->GetNumPoints();
299  const int nquad2 = m_base[2]->GetNumPoints();
300  const int order0 = m_base[0]->GetNumModes();
301  const int order1 = m_base[1]->GetNumModes();
302  Array<OneD, NekDouble> wsp(nquad1 * nquad2 * order0 +
303  nquad2 * order0 * (order1 + 1) / 2);
304 
305  if (multiplybyweights)
306  {
307  Array<OneD, NekDouble> tmp(nquad0 * nquad1 * nquad2);
308 
309  MultiplyByQuadratureMetric(inarray, tmp);
311  m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
312  tmp, outarray, wsp, true, true, true);
313  }
314  else
315  {
317  m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
318  inarray, outarray, wsp, true, true, true);
319  }
320 }
void IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)

References Nektar::StdRegions::StdExpansion3D::IProductWRTBase_SumFacKernel(), Nektar::StdRegions::StdExpansion::m_base, and Nektar::StdRegions::StdExpansion::MultiplyByQuadratureMetric().

Referenced by v_IProductWRTBase().

◆ v_IProductWRTDerivBase()

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

Calculates the inner product \( I_{pqr} = (u, \partial_{x_i} \phi_{pqr}) \).

The derivative of the basis functions is performed using the chain rule in order to incorporate the geometric factors. Assuming that the basis functions are a tensor product \(\phi_{pqr}(\eta_1,\eta_2,\eta_3) = \phi_1(\eta_1)\phi_2(\eta_2)\phi_3(\eta_3)\), this yields the result

\[ I_{pqr} = \sum_{j=1}^3 \left(u, \frac{\partial u}{\partial \eta_j} \frac{\partial \eta_j}{\partial x_i}\right) \]

In the prismatic element, we must also incorporate a second set of geometric factors which incorporate the collapsed co-ordinate system, so that

\[ \frac{\partial\eta_j}{\partial x_i} = \sum_{k=1}^3 \frac{\partial\eta_j}{\partial\xi_k}\frac{\partial\xi_k}{\partial x_i} \]

These derivatives can be found on p152 of Sherwin & Karniadakis.

Parameters
dirDirection in which to take the derivative.
inarrayThe function \( u \).
outarrayValue of the inner product.

Reimplemented from Nektar::StdRegions::StdTetExp.

Definition at line 352 of file TetExp.cpp.

355 {
356  const int nquad0 = m_base[0]->GetNumPoints();
357  const int nquad1 = m_base[1]->GetNumPoints();
358  const int nquad2 = m_base[2]->GetNumPoints();
359  const int order0 = m_base[0]->GetNumModes();
360  const int order1 = m_base[1]->GetNumModes();
361  const int nqtot = nquad0 * nquad1 * nquad2;
362 
363  Array<OneD, NekDouble> tmp1(nqtot);
364  Array<OneD, NekDouble> tmp2(nqtot);
365  Array<OneD, NekDouble> tmp3(nqtot);
366  Array<OneD, NekDouble> tmp4(nqtot);
367  Array<OneD, NekDouble> tmp6(m_ncoeffs);
368  Array<OneD, NekDouble> wsp(nquad1 * nquad2 * order0 +
369  nquad2 * order0 * (order1 + 1) / 2);
370 
371  MultiplyByQuadratureMetric(inarray, tmp1);
372 
373  Array<OneD, Array<OneD, NekDouble>> tmp2D{3};
374  tmp2D[0] = tmp2;
375  tmp2D[1] = tmp3;
376  tmp2D[2] = tmp4;
377 
378  TetExp::v_AlignVectorToCollapsedDir(dir, tmp1, tmp2D);
379 
380  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
381  m_base[2]->GetBdata(), tmp2, outarray, wsp,
382  false, true, true);
383 
384  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
385  m_base[2]->GetBdata(), tmp3, tmp6, wsp, true,
386  false, true);
387 
388  Vmath::Vadd(m_ncoeffs, tmp6, 1, outarray, 1, outarray, 1);
389 
390  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetBdata(),
391  m_base[2]->GetDbdata(), tmp4, tmp6, wsp, true,
392  true, false);
393 
394  Vmath::Vadd(m_ncoeffs, tmp6, 1, outarray, 1, outarray, 1);
395 }
virtual void v_AlignVectorToCollapsedDir(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
Definition: TetExp.cpp:397

References Nektar::StdRegions::StdExpansion3D::IProductWRTBase_SumFacKernel(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::StdRegions::StdExpansion::MultiplyByQuadratureMetric(), v_AlignVectorToCollapsedDir(), and Vmath::Vadd().

◆ v_LaplacianMatrixOp() [1/2]

void Nektar::LocalRegions::TetExp::v_LaplacianMatrixOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdRegions::StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 975 of file TetExp.cpp.

978 {
979  TetExp::v_LaplacianMatrixOp_MatFree(inarray, outarray, mkey);
980 }
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)

References Nektar::StdRegions::StdExpansion3D::v_LaplacianMatrixOp_MatFree().

◆ v_LaplacianMatrixOp() [2/2]

void Nektar::LocalRegions::TetExp::v_LaplacianMatrixOp ( const int  k1,
const int  k2,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdRegions::StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 982 of file TetExp.cpp.

986 {
987  StdExpansion::LaplacianMatrixOp_MatFree(k1, k2, inarray, outarray, mkey);
988 }

◆ v_LaplacianMatrixOp_MatFree_Kernel()

void Nektar::LocalRegions::TetExp::v_LaplacianMatrixOp_MatFree_Kernel ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
Array< OneD, NekDouble > &  wsp 
)
privatevirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1091 of file TetExp.cpp.

1094 {
1095  // This implementation is only valid when there are no
1096  // coefficients associated to the Laplacian operator
1097  if (m_metrics.count(eMetricLaplacian00) == 0)
1098  {
1100  }
1101 
1102  int nquad0 = m_base[0]->GetNumPoints();
1103  int nquad1 = m_base[1]->GetNumPoints();
1104  int nquad2 = m_base[2]->GetNumPoints();
1105  int nqtot = nquad0 * nquad1 * nquad2;
1106 
1107  ASSERTL1(wsp.size() >= 6 * nqtot, "Insufficient workspace size.");
1108  ASSERTL1(m_ncoeffs <= nqtot, "Workspace not set up for ncoeffs > nqtot");
1109 
1110  const Array<OneD, const NekDouble> &base0 = m_base[0]->GetBdata();
1111  const Array<OneD, const NekDouble> &base1 = m_base[1]->GetBdata();
1112  const Array<OneD, const NekDouble> &base2 = m_base[2]->GetBdata();
1113  const Array<OneD, const NekDouble> &dbase0 = m_base[0]->GetDbdata();
1114  const Array<OneD, const NekDouble> &dbase1 = m_base[1]->GetDbdata();
1115  const Array<OneD, const NekDouble> &dbase2 = m_base[2]->GetDbdata();
1116  const Array<OneD, const NekDouble> &metric00 =
1117  m_metrics[eMetricLaplacian00];
1118  const Array<OneD, const NekDouble> &metric01 =
1119  m_metrics[eMetricLaplacian01];
1120  const Array<OneD, const NekDouble> &metric02 =
1121  m_metrics[eMetricLaplacian02];
1122  const Array<OneD, const NekDouble> &metric11 =
1123  m_metrics[eMetricLaplacian11];
1124  const Array<OneD, const NekDouble> &metric12 =
1125  m_metrics[eMetricLaplacian12];
1126  const Array<OneD, const NekDouble> &metric22 =
1127  m_metrics[eMetricLaplacian22];
1128 
1129  // Allocate temporary storage
1130  Array<OneD, NekDouble> wsp0(2 * nqtot, wsp);
1131  Array<OneD, NekDouble> wsp1(nqtot, wsp + 1 * nqtot);
1132  Array<OneD, NekDouble> wsp2(nqtot, wsp + 2 * nqtot);
1133  Array<OneD, NekDouble> wsp3(nqtot, wsp + 3 * nqtot);
1134  Array<OneD, NekDouble> wsp4(nqtot, wsp + 4 * nqtot);
1135  Array<OneD, NekDouble> wsp5(nqtot, wsp + 5 * nqtot);
1136 
1137  // LAPLACIAN MATRIX OPERATION
1138  // wsp1 = du_dxi1 = D_xi1 * inarray = D_xi1 * u
1139  // wsp2 = du_dxi2 = D_xi2 * inarray = D_xi2 * u
1140  // wsp2 = du_dxi3 = D_xi3 * inarray = D_xi3 * u
1141  StdExpansion3D::PhysTensorDeriv(inarray, wsp0, wsp1, wsp2);
1142 
1143  // wsp0 = k = g0 * wsp1 + g1 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
1144  // wsp2 = l = g1 * wsp1 + g2 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
1145  // where g0, g1 and g2 are the metric terms set up in the GeomFactors class
1146  // especially for this purpose
1147  Vmath::Vvtvvtp(nqtot, &metric00[0], 1, &wsp0[0], 1, &metric01[0], 1,
1148  &wsp1[0], 1, &wsp3[0], 1);
1149  Vmath::Vvtvp(nqtot, &metric02[0], 1, &wsp2[0], 1, &wsp3[0], 1, &wsp3[0], 1);
1150  Vmath::Vvtvvtp(nqtot, &metric01[0], 1, &wsp0[0], 1, &metric11[0], 1,
1151  &wsp1[0], 1, &wsp4[0], 1);
1152  Vmath::Vvtvp(nqtot, &metric12[0], 1, &wsp2[0], 1, &wsp4[0], 1, &wsp4[0], 1);
1153  Vmath::Vvtvvtp(nqtot, &metric02[0], 1, &wsp0[0], 1, &metric12[0], 1,
1154  &wsp1[0], 1, &wsp5[0], 1);
1155  Vmath::Vvtvp(nqtot, &metric22[0], 1, &wsp2[0], 1, &wsp5[0], 1, &wsp5[0], 1);
1156 
1157  // outarray = m = (D_xi1 * B)^T * k
1158  // wsp1 = n = (D_xi2 * B)^T * l
1159  IProductWRTBase_SumFacKernel(dbase0, base1, base2, wsp3, outarray, wsp0,
1160  false, true, true);
1161  IProductWRTBase_SumFacKernel(base0, dbase1, base2, wsp4, wsp2, wsp0, true,
1162  false, true);
1163  Vmath::Vadd(m_ncoeffs, wsp2.get(), 1, outarray.get(), 1, outarray.get(), 1);
1164  IProductWRTBase_SumFacKernel(base0, base1, dbase2, wsp5, wsp2, wsp0, true,
1165  true, false);
1166  Vmath::Vadd(m_ncoeffs, wsp2.get(), 1, outarray.get(), 1, outarray.get(), 1);
1167 }

References ASSERTL1, Nektar::LocalRegions::Expansion::ComputeLaplacianMetric(), Nektar::LocalRegions::eMetricLaplacian00, Nektar::StdRegions::StdExpansion::m_base, and Nektar::LocalRegions::Expansion::m_metrics.

◆ v_NormalTraceDerivFactors()

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

: This method gets all of the factors which are required as part of the Gradient Jump Penalty stabilisation and involves the product of the normal and geometric factors along the element trace.

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1372 of file TetExp.cpp.

1376 {
1377  int nquad0 = GetNumPoints(0);
1378  int nquad1 = GetNumPoints(1);
1379  int nquad2 = GetNumPoints(2);
1380 
1381  const Array<TwoD, const NekDouble> &df =
1382  m_metricinfo->GetDerivFactors(GetPointsKeys());
1383 
1384  if (d0factors.size() != 4)
1385  {
1386  d0factors = Array<OneD, Array<OneD, NekDouble>>(4);
1387  d1factors = Array<OneD, Array<OneD, NekDouble>>(4);
1388  d2factors = Array<OneD, Array<OneD, NekDouble>>(4);
1389  }
1390 
1391  if (d0factors[0].size() != nquad0 * nquad1)
1392  {
1393  d0factors[0] = Array<OneD, NekDouble>(nquad0 * nquad1);
1394  d1factors[0] = Array<OneD, NekDouble>(nquad0 * nquad1);
1395  d2factors[0] = Array<OneD, NekDouble>(nquad0 * nquad1);
1396  }
1397 
1398  if (d0factors[1].size() != nquad0 * nquad2)
1399  {
1400  d0factors[1] = Array<OneD, NekDouble>(nquad0 * nquad2);
1401  d1factors[1] = Array<OneD, NekDouble>(nquad0 * nquad2);
1402  d2factors[1] = Array<OneD, NekDouble>(nquad0 * nquad2);
1403  }
1404 
1405  if (d0factors[2].size() != nquad1 * nquad2)
1406  {
1407  d0factors[2] = Array<OneD, NekDouble>(nquad1 * nquad2);
1408  d0factors[3] = Array<OneD, NekDouble>(nquad1 * nquad2);
1409  d1factors[2] = Array<OneD, NekDouble>(nquad1 * nquad2);
1410  d1factors[3] = Array<OneD, NekDouble>(nquad1 * nquad2);
1411  d2factors[2] = Array<OneD, NekDouble>(nquad1 * nquad2);
1412  d2factors[3] = Array<OneD, NekDouble>(nquad1 * nquad2);
1413  }
1414 
1415  // Outwards normals
1416  const Array<OneD, const Array<OneD, NekDouble>> &normal_0 =
1417  GetTraceNormal(0);
1418  const Array<OneD, const Array<OneD, NekDouble>> &normal_1 =
1419  GetTraceNormal(1);
1420  const Array<OneD, const Array<OneD, NekDouble>> &normal_2 =
1421  GetTraceNormal(2);
1422  const Array<OneD, const Array<OneD, NekDouble>> &normal_3 =
1423  GetTraceNormal(3);
1424 
1425  int ncoords = normal_0.size();
1426 
1427  // first gather together standard cartesian inner products
1428  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1429  {
1430  // face 0
1431  for (int i = 0; i < nquad0 * nquad1; ++i)
1432  {
1433  d0factors[0][i] = df[0][i] * normal_0[0][i];
1434  d1factors[0][i] = df[1][i] * normal_0[0][i];
1435  d2factors[0][i] = df[2][i] * normal_0[0][i];
1436  }
1437 
1438  for (int n = 1; n < ncoords; ++n)
1439  {
1440  for (int i = 0; i < nquad0 * nquad1; ++i)
1441  {
1442  d0factors[0][i] += df[3 * n][i] * normal_0[n][i];
1443  d1factors[0][i] += df[3 * n + 1][i] * normal_0[n][i];
1444  d2factors[0][i] += df[3 * n + 2][i] * normal_0[n][i];
1445  }
1446  }
1447 
1448  // face 1
1449  for (int j = 0; j < nquad2; ++j)
1450  {
1451  for (int i = 0; i < nquad0; ++i)
1452  {
1453  d0factors[1][i] = df[0][j * nquad0 * nquad1 + i] *
1454  normal_1[0][j * nquad0 + i];
1455  d1factors[1][i] = df[1][j * nquad0 * nquad1 + i] *
1456  normal_1[0][j * nquad0 + i];
1457  d2factors[1][i] = df[2][j * nquad0 * nquad1 + i] *
1458  normal_1[0][j * nquad0 + i];
1459  }
1460  }
1461 
1462  for (int n = 1; n < ncoords; ++n)
1463  {
1464  for (int j = 0; j < nquad2; ++j)
1465  {
1466  for (int i = 0; i < nquad0; ++i)
1467  {
1468  d0factors[1][i] = df[3 * n][j * nquad0 * nquad1 + i] *
1469  normal_1[0][j * nquad0 + i];
1470  d1factors[1][i] = df[3 * n + 1][j * nquad0 * nquad1 + i] *
1471  normal_1[0][j * nquad0 + i];
1472  d2factors[1][i] = df[3 * n + 2][j * nquad0 * nquad1 + i] *
1473  normal_1[0][j * nquad0 + i];
1474  }
1475  }
1476  }
1477 
1478  // faces 2 and 3
1479  for (int j = 0; j < nquad2; ++j)
1480  {
1481  for (int i = 0; i < nquad1; ++i)
1482  {
1483  d0factors[2][j * nquad1 + i] =
1484  df[0][j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1485  normal_2[0][j * nquad1 + i];
1486  d0factors[3][j * nquad1 + i] =
1487  df[0][j * nquad0 * nquad1 + i * nquad0] *
1488  normal_3[0][j * nquad1 + i];
1489  d1factors[2][j * nquad1 + i] =
1490  df[1][j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1491  normal_2[0][j * nquad1 + i];
1492  d1factors[3][j * nquad1 + i] =
1493  df[1][j * nquad0 * nquad1 + i * nquad0] *
1494  normal_3[0][j * nquad1 + i];
1495  d2factors[2][j * nquad1 + i] =
1496  df[2][j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1497  normal_2[0][j * nquad1 + i];
1498  d2factors[3][j * nquad1 + i] =
1499  df[2][j * nquad0 * nquad1 + i * nquad0] *
1500  normal_3[0][j * nquad1 + i];
1501  }
1502  }
1503 
1504  for (int n = 1; n < ncoords; ++n)
1505  {
1506  for (int j = 0; j < nquad2; ++j)
1507  {
1508  for (int i = 0; i < nquad1; ++i)
1509  {
1510  d0factors[2][j * nquad1 + i] +=
1511  df[3 * n][j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1512  normal_2[n][j * nquad0 + i];
1513  d0factors[3][j * nquad0 + i] +=
1514  df[3 * n][i * nquad0 + j * nquad0 * nquad1] *
1515  normal_3[n][j * nquad0 + i];
1516  d1factors[2][j * nquad1 + i] +=
1517  df[3 * n + 1]
1518  [j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1519  normal_2[n][j * nquad0 + i];
1520  d1factors[3][j * nquad0 + i] +=
1521  df[3 * n + 1][i * nquad0 + j * nquad0 * nquad1] *
1522  normal_3[n][j * nquad0 + i];
1523  d2factors[2][j * nquad1 + i] +=
1524  df[3 * n + 2]
1525  [j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1526  normal_2[n][j * nquad0 + i];
1527  d2factors[3][j * nquad0 + i] +=
1528  df[3 * n + 2][i * nquad0 + j * nquad0 * nquad1] *
1529  normal_3[n][j * nquad0 + i];
1530  }
1531  }
1532  }
1533  }
1534  else
1535  {
1536  // Face 0
1537  for (int i = 0; i < nquad0 * nquad1; ++i)
1538  {
1539  d0factors[0][i] = df[0][0] * normal_0[0][i];
1540  d1factors[0][i] = df[1][0] * normal_0[0][i];
1541  d2factors[0][i] = df[2][0] * normal_0[0][i];
1542  }
1543 
1544  for (int n = 1; n < ncoords; ++n)
1545  {
1546  for (int i = 0; i < nquad0 * nquad1; ++i)
1547  {
1548  d0factors[0][i] += df[3 * n][0] * normal_0[n][i];
1549  d1factors[0][i] += df[3 * n + 1][0] * normal_0[n][i];
1550  d2factors[0][i] += df[3 * n + 2][0] * normal_0[n][i];
1551  }
1552  }
1553 
1554  // face 1
1555  for (int i = 0; i < nquad0 * nquad2; ++i)
1556  {
1557  d0factors[1][i] = df[0][0] * normal_1[0][i];
1558  d1factors[1][i] = df[1][0] * normal_1[0][i];
1559  d2factors[1][i] = df[2][0] * normal_1[0][i];
1560  }
1561 
1562  for (int n = 1; n < ncoords; ++n)
1563  {
1564  for (int i = 0; i < nquad0 * nquad2; ++i)
1565  {
1566  d0factors[1][i] += df[3 * n][0] * normal_1[n][i];
1567  d1factors[1][i] += df[3 * n + 1][0] * normal_1[n][i];
1568  d2factors[1][i] += df[3 * n + 2][0] * normal_1[n][i];
1569  }
1570  }
1571 
1572  // faces 2 and 3
1573  for (int i = 0; i < nquad1 * nquad2; ++i)
1574  {
1575  d0factors[2][i] = df[0][0] * normal_2[0][i];
1576  d0factors[3][i] = df[0][0] * normal_3[0][i];
1577 
1578  d1factors[2][i] = df[1][0] * normal_2[0][i];
1579  d1factors[3][i] = df[1][0] * normal_3[0][i];
1580 
1581  d2factors[2][i] = df[2][0] * normal_2[0][i];
1582  d2factors[3][i] = df[2][0] * normal_3[0][i];
1583  }
1584 
1585  for (int n = 1; n < ncoords; ++n)
1586  {
1587  for (int i = 0; i < nquad1 * nquad2; ++i)
1588  {
1589  d0factors[2][i] += df[3 * n][0] * normal_2[n][i];
1590  d0factors[3][i] += df[3 * n][0] * normal_3[n][i];
1591 
1592  d1factors[2][i] += df[3 * n + 1][0] * normal_2[n][i];
1593  d1factors[3][i] += df[3 * n + 1][0] * normal_3[n][i];
1594 
1595  d2factors[2][i] += df[3 * n + 2][0] * normal_2[n][i];
1596  d2factors[3][i] += df[3 * n + 2][0] * normal_3[n][i];
1597  }
1598  }
1599  }
1600 }
const NormalVector & GetTraceNormal(const int id)
Definition: Expansion.cpp:250

◆ v_PhysDeriv()

void Nektar::LocalRegions::TetExp::v_PhysDeriv ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d0,
Array< OneD, NekDouble > &  out_d1,
Array< OneD, NekDouble > &  out_d2 
)
protectedvirtual

Differentiate inarray in the three coordinate directions.

Parameters
inarrayInput array of values at quadrature points to be differentiated.
out_d0Derivative in first coordinate direction.
out_d1Derivative in second coordinate direction.
out_d2Derivative in third coordinate direction.

Reimplemented from Nektar::StdRegions::StdTetExp.

Definition at line 152 of file TetExp.cpp.

156 {
157  int TotPts = m_base[0]->GetNumPoints() * m_base[1]->GetNumPoints() *
158  m_base[2]->GetNumPoints();
159 
160  Array<TwoD, const NekDouble> df =
161  m_metricinfo->GetDerivFactors(GetPointsKeys());
162  Array<OneD, NekDouble> Diff0 = Array<OneD, NekDouble>(3 * TotPts);
163  Array<OneD, NekDouble> Diff1 = Diff0 + TotPts;
164  Array<OneD, NekDouble> Diff2 = Diff1 + TotPts;
165 
166  StdTetExp::v_PhysDeriv(inarray, Diff0, Diff1, Diff2);
167 
168  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
169  {
170  if (out_d0.size())
171  {
172  Vmath::Vmul(TotPts, &df[0][0], 1, &Diff0[0], 1, &out_d0[0], 1);
173  Vmath::Vvtvp(TotPts, &df[1][0], 1, &Diff1[0], 1, &out_d0[0], 1,
174  &out_d0[0], 1);
175  Vmath::Vvtvp(TotPts, &df[2][0], 1, &Diff2[0], 1, &out_d0[0], 1,
176  &out_d0[0], 1);
177  }
178 
179  if (out_d1.size())
180  {
181  Vmath::Vmul(TotPts, &df[3][0], 1, &Diff0[0], 1, &out_d1[0], 1);
182  Vmath::Vvtvp(TotPts, &df[4][0], 1, &Diff1[0], 1, &out_d1[0], 1,
183  &out_d1[0], 1);
184  Vmath::Vvtvp(TotPts, &df[5][0], 1, &Diff2[0], 1, &out_d1[0], 1,
185  &out_d1[0], 1);
186  }
187 
188  if (out_d2.size())
189  {
190  Vmath::Vmul(TotPts, &df[6][0], 1, &Diff0[0], 1, &out_d2[0], 1);
191  Vmath::Vvtvp(TotPts, &df[7][0], 1, &Diff1[0], 1, &out_d2[0], 1,
192  &out_d2[0], 1);
193  Vmath::Vvtvp(TotPts, &df[8][0], 1, &Diff2[0], 1, &out_d2[0], 1,
194  &out_d2[0], 1);
195  }
196  }
197  else // regular geometry
198  {
199  if (out_d0.size())
200  {
201  Vmath::Smul(TotPts, df[0][0], &Diff0[0], 1, &out_d0[0], 1);
202  Blas::Daxpy(TotPts, df[1][0], &Diff1[0], 1, &out_d0[0], 1);
203  Blas::Daxpy(TotPts, df[2][0], &Diff2[0], 1, &out_d0[0], 1);
204  }
205 
206  if (out_d1.size())
207  {
208  Vmath::Smul(TotPts, df[3][0], &Diff0[0], 1, &out_d1[0], 1);
209  Blas::Daxpy(TotPts, df[4][0], &Diff1[0], 1, &out_d1[0], 1);
210  Blas::Daxpy(TotPts, df[5][0], &Diff2[0], 1, &out_d1[0], 1);
211  }
212 
213  if (out_d2.size())
214  {
215  Vmath::Smul(TotPts, df[6][0], &Diff0[0], 1, &out_d2[0], 1);
216  Blas::Daxpy(TotPts, df[7][0], &Diff1[0], 1, &out_d2[0], 1);
217  Blas::Daxpy(TotPts, df[8][0], &Diff2[0], 1, &out_d2[0], 1);
218  }
219  }
220 }
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
Definition: Blas.hpp:154

References Blas::Daxpy(), Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_metricinfo, Vmath::Smul(), Vmath::Vmul(), and Vmath::Vvtvp().

◆ v_PhysEvaluate()

NekDouble Nektar::LocalRegions::TetExp::v_PhysEvaluate ( const Array< OneD, const NekDouble > &  coord,
const Array< OneD, const NekDouble > &  physvals 
)
protectedvirtual
Parameters
coordPhysical space coordinate
Returns
Evaluation of expansion at given coordinate.

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 480 of file TetExp.cpp.

482 {
483  ASSERTL0(m_geom, "m_geom not defined");
484 
485  Array<OneD, NekDouble> Lcoord = Array<OneD, NekDouble>(3);
486 
487  // Get the local (eta) coordinates of the point
488  m_geom->GetLocCoords(coord, Lcoord);
489 
490  // Evaluate point in local (eta) coordinates.
491  return StdTetExp::v_PhysEvaluate(Lcoord, physvals);
492 }

References ASSERTL0, and Nektar::LocalRegions::Expansion::m_geom.

◆ v_StdPhysEvaluate()

NekDouble Nektar::LocalRegions::TetExp::v_StdPhysEvaluate ( const Array< OneD, const NekDouble > &  Lcoord,
const Array< OneD, const NekDouble > &  physvals 
)
protectedvirtual

Given the local cartesian coordinate Lcoord evaluate the value of physvals at this point by calling through to the StdExpansion method

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 468 of file TetExp.cpp.

471 {
472  // Evaluate point in local (eta) coordinates.
473  return StdTetExp::v_PhysEvaluate(Lcoord, physvals);
474 }

◆ v_SVVLaplacianFilter()

void Nektar::LocalRegions::TetExp::v_SVVLaplacianFilter ( Array< OneD, NekDouble > &  array,
const StdRegions::StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdTetExp.

Definition at line 990 of file TetExp.cpp.

992 {
993  int nq = GetTotPoints();
994 
995  // Calculate sqrt of the Jacobian
996  Array<OneD, const NekDouble> jac = m_metricinfo->GetJac(GetPointsKeys());
997  Array<OneD, NekDouble> sqrt_jac(nq);
998  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
999  {
1000  Vmath::Vsqrt(nq, jac, 1, sqrt_jac, 1);
1001  }
1002  else
1003  {
1004  Vmath::Fill(nq, sqrt(jac[0]), sqrt_jac, 1);
1005  }
1006 
1007  // Multiply array by sqrt(Jac)
1008  Vmath::Vmul(nq, sqrt_jac, 1, array, 1, array, 1);
1009 
1010  // Apply std region filter
1011  StdTetExp::v_SVVLaplacianFilter(array, mkey);
1012 
1013  // Divide by sqrt(Jac)
1014  Vmath::Vdiv(nq, array, 1, sqrt_jac, 1, array, 1);
1015 }
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 Nektar::SpatialDomains::eDeformed, Vmath::Fill(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::GetTotPoints(), Nektar::LocalRegions::Expansion::m_metricinfo, tinysimd::sqrt(), Vmath::Vdiv(), Vmath::Vmul(), and Vmath::Vsqrt().

Member Data Documentation

◆ m_matrixManager

LibUtilities::NekManager<MatrixKey, DNekScalMat, MatrixKey::opLess> Nektar::LocalRegions::TetExp::m_matrixManager
private

Definition at line 194 of file TetExp.h.

Referenced by v_FwdTrans(), and v_GetLocMatrix().

◆ m_staticCondMatrixManager

LibUtilities::NekManager<MatrixKey, DNekScalBlkMat, MatrixKey::opLess> Nektar::LocalRegions::TetExp::m_staticCondMatrixManager
private

Definition at line 196 of file TetExp.h.

Referenced by v_DropLocStaticCondMatrix(), and v_GetLocStaticCondMatrix().