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

#include <HexExp.h>

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

Public Member Functions

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

Protected Member Functions

NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray) override
 Integrate the physical point list inarray over region. More...
 
void v_PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2) override
 Calculate the derivative of the physical points. More...
 
void v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Calculate the derivative of the physical points in a single direction. More...
 
void v_PhysDirectionalDeriv (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &out) override
 Physical derivative along a direction vector. More...
 
void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Forward transform from physical quadrature space stored in inarray and evaluate the expansion coefficients and store in (this)->_coeffs. More...
 
void v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Calculate the inner product of inarray with respect to the elements basis. More...
 
void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
 Calculate the inner product of inarray with respect to the given basis B = base0 * base1 * base2. More...
 
void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 Calculates the inner product \( I_{pqr} = (u, \partial_{x_i} \phi_{pqr}) \). More...
 
void v_AlignVectorToCollapsedDir (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
 
void v_IProductWRTDirectionalDerivBase (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_IProductWRTDirectionalDerivBase_SumFac (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals) override
 
NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals) override
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble v_PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
 
void v_GetCoord (const Array< OneD, const NekDouble > &Lcoords, Array< OneD, NekDouble > &coords) override
 Retrieves the physical coordinates of a given set of reference coordinates. More...
 
void v_GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
 
LibUtilities::ShapeType v_DetShapeType () const override
 Return the region shape using the enum-list of ShapeType. More...
 
StdRegions::StdExpansionSharedPtr v_GetStdExp (void) const override
 
StdRegions::StdExpansionSharedPtr v_GetLinStdExp (void) const override
 
void v_ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int mode_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType) override
 
void v_GetTracePhysMap (const int face, Array< OneD, int > &outarray) override
 
void v_ComputeTraceNormal (const int face) override
 
void v_MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
void v_LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
void v_LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
void v_WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
void v_WeakDirectionalDerivMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
void v_MassLevelCurvatureMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
void v_HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
void v_ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdRegions::StdMatrixKey &mkey) override
 
DNekMatSharedPtr v_GenMatrix (const StdRegions::StdMatrixKey &mkey) override
 
DNekMatSharedPtr v_CreateStdMatrix (const StdRegions::StdMatrixKey &mkey) override
 
DNekScalMatSharedPtr v_GetLocMatrix (const MatrixKey &mkey) override
 
void v_DropLocMatrix (const MatrixKey &mkey) override
 
DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const MatrixKey &mkey) override
 
void v_DropLocStaticCondMatrix (const MatrixKey &mkey) override
 
void v_ComputeLaplacianMetric () override
 
- Protected Member Functions inherited from Nektar::StdRegions::StdHexExp
void v_PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2) override
 Differentiation Methods. More...
 
void v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2) override
 
void v_StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
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) override
 
void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multbyweights=true) override
 
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) override
 
void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta) override
 
void v_LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi) override
 
void v_FillMode (const int mode, Array< OneD, NekDouble > &outarray) override
 
NekDouble v_PhysEvaluateBasis (const Array< OneD, const NekDouble > &coords, int mode) final
 
NekDouble v_PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
 
int v_GetNverts () const override
 
int v_GetNedges () const override
 
int v_GetNtraces () const override
 
LibUtilities::ShapeType v_DetShapeType () const override
 
int v_NumBndryCoeffs () const override
 
int v_NumDGBndryCoeffs () const override
 
int v_GetTraceNcoeffs (const int i) const override
 
int v_GetTraceIntNcoeffs (const int i) const override
 
int v_GetTraceNumPoints (const int i) const override
 
LibUtilities::PointsKey v_GetTracePointsKey (const int i, const int j) const override
 
int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset) override
 
const LibUtilities::BasisKey v_GetTraceBasisKey (const int i, const int k) const override
 
bool v_IsBoundaryInteriorExpansion () const override
 
void v_GetCoords (Array< OneD, NekDouble > &coords_x, Array< OneD, NekDouble > &coords_y, Array< OneD, NekDouble > &coords_z) override
 
void v_GetTraceNumModes (const int fid, int &numModes0, int &numModes1, Orientation faceOrient=eDir1FwdDir1_Dir2FwdDir2) override
 
int v_GetEdgeNcoeffs (const int i) const override
 
int v_GetVertexMap (int localVertexId, bool useCoeffPacking=false) override
 
void v_GetInteriorMap (Array< OneD, unsigned int > &outarray) override
 
void v_GetBoundaryMap (Array< OneD, unsigned int > &outarray) override
 
void v_GetTraceCoeffMap (const unsigned int fid, Array< OneD, unsigned int > &maparray) override
 
void v_GetElmtTraceToTraceMap (const unsigned int fid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation faceOrient, int P, int Q) override
 
void v_GetEdgeInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2) override
 
void v_GetTraceInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2) override
 
DNekMatSharedPtr v_GenMatrix (const StdMatrixKey &mkey) override
 
DNekMatSharedPtr v_CreateStdMatrix (const StdMatrixKey &mkey) override
 
void v_GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true) override
 
void v_MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
 
void v_LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
 
void v_LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
 
void v_WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
 
void v_HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey) override
 
void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey) override
 
void v_ExponentialFilter (Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff) override
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion3D
NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals) override
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble v_PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals) override
 
NekDouble v_PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
 
virtual void v_BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)=0
 
virtual void v_IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)=0
 
void v_LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
void v_HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray) override
 Integrates the specified function over the domain. More...
 
virtual int v_GetNedges (void) const
 
virtual int v_GetEdgeNcoeffs (const int i) const
 
NekDouble BaryTensorDeriv (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
 
virtual void v_GetEdgeInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards)
 
void v_GetTraceToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient, int P, int Q) override
 
void v_GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat) override
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion
DNekMatSharedPtr CreateStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr CreateStdStaticCondMatrix (const StdMatrixKey &mkey)
 Create the static condensation of a matrix when using a boundary interior decomposition. More...
 
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase_SumFac (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GeneralMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
 
void LaplacianMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp_MatFree (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void HelmholtzMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
 
virtual NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
virtual void v_GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 
virtual void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals, NekDouble &deriv, NekDouble &deriv2)
 This function performs the barycentric interpolation of the polynomial stored in coord at a point physvals using barycentric interpolation weights in direction. More...
 
template<int DIR>
NekDouble BaryEvaluateBasis (const NekDouble &coord, const int &mode)
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals)
 Helper function to pass an unused value by reference into BaryEvaluate. More...
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals, NekDouble &deriv)
 
- Protected Member Functions inherited from Nektar::LocalRegions::Expansion3D
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) override
 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...
 
DNekMatSharedPtr v_GenMatrix (const StdRegions::StdMatrixKey &mkey) override
 
void v_AddFaceNormBoundaryInt (const int face, const ExpansionSharedPtr &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray) override
 
void v_AddRobinMassMatrix (const int face, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat) override
 
StdRegions::Orientation v_GetTraceOrient (int face) override
 
void v_GetTracePhysVals (const int face, const StdRegions::StdExpansionSharedPtr &FaceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient) override
 Extract the physical values along face face from inarray into outarray following the local face orientation and point distribution defined by defined in FaceExp. More...
 
void v_GenTraceExp (const int traceid, ExpansionSharedPtr &exp) override
 
void GetPhysFaceVarCoeffsFromElement (const int face, ExpansionSharedPtr &FaceExp, const Array< OneD, const NekDouble > &varcoeff, Array< OneD, NekDouble > &outarray)
 
DNekMatSharedPtr v_BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType) override
 
DNekMatSharedPtr v_BuildInverseTransformationMatrix (const DNekScalMatSharedPtr &transformationmatrix) override
 Build inverse and inverse transposed transformation matrix: \(\mathbf{R^{-1}}\) and \(\mathbf{R^{-T}}\). More...
 
DNekMatSharedPtr v_BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd) override
 
void v_TraceNormLen (const int traceid, NekDouble &h, NekDouble &p) override
 
- 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)
 
Array< OneD, NekDoubleGetMF (const int dir, const int shapedim, const StdRegions::VarCoeffMap &varcoeffs)
 
Array< OneD, NekDoubleGetMFDiv (const int dir, const StdRegions::VarCoeffMap &varcoeffs)
 
Array< OneD, NekDoubleGetMFMag (const int dir, const StdRegions::VarCoeffMap &varcoeffs)
 
void v_MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
virtual void v_DivideByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_ComputeLaplacianMetric ()
 
int v_GetCoordim () const override
 
void v_GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
 
virtual DNekScalMatSharedPtr v_GetLocMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual void v_DropLocMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual DNekMatSharedPtr v_BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
 
virtual DNekMatSharedPtr v_BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd)
 
virtual void v_ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType)
 
virtual void v_AddEdgeNormBoundaryInt (const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddEdgeNormBoundaryInt (const int edge, const std::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddFaceNormBoundaryInt (const int face, const std::shared_ptr< Expansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
virtual void v_DGDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &coeffs, Array< OneD, NekDouble > &outarray)
 
virtual NekDouble v_VectorFlux (const Array< OneD, Array< OneD, NekDouble > > &vec)
 
virtual void v_NormalTraceDerivFactors (Array< OneD, Array< OneD, NekDouble > > &factors, Array< OneD, Array< OneD, NekDouble > > &d0factors, Array< OneD, Array< OneD, NekDouble > > &d1factors)
 
virtual void v_AlignVectorToCollapsedDir (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual StdRegions::Orientation v_GetTraceOrient (int trace)
 
void v_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
virtual void v_GetTraceQFactors (const int trace, Array< OneD, NekDouble > &outarray)
 
virtual void v_GetTracePhysVals (const int trace, const StdRegions::StdExpansionSharedPtr &TraceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient)
 
virtual void v_GetTracePhysMap (const int edge, Array< OneD, int > &outarray)
 
virtual void v_ReOrientTracePhysMap (const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1=-1)
 
virtual void v_ComputeTraceNormal (const int id)
 
virtual const Array< OneD, const NekDouble > & v_GetPhysNormals ()
 
virtual void v_SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
virtual void v_SetUpPhysNormals (const int id)
 
virtual void v_AddRobinMassMatrix (const int face, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
virtual void v_AddRobinTraceContribution (const int traceid, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs)
 
virtual void v_TraceNormLen (const int traceid, NekDouble &h, NekDouble &p)
 
virtual void v_GenTraceExp (const int traceid, ExpansionSharedPtr &exp)
 

Private Member Functions

void v_LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp) override
 
void v_NormalTraceDerivFactors (Array< OneD, Array< OneD, NekDouble > > &factors, Array< OneD, Array< OneD, NekDouble > > &d0factors, Array< OneD, Array< OneD, NekDouble > > &d1factors) override
 : 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...
 

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 hexahedral local expansion.

Definition at line 48 of file HexExp.h.

Constructor & Destructor Documentation

◆ HexExp() [1/2]

Nektar::LocalRegions::HexExp::HexExp ( const LibUtilities::BasisKey Ba,
const LibUtilities::BasisKey Bb,
const LibUtilities::BasisKey Bc,
const SpatialDomains::HexGeomSharedPtr 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 57 of file HexExp.cpp.

61 : StdExpansion(Ba.GetNumModes() * Bb.GetNumModes() * Bc.GetNumModes(), 3,
62 Ba, Bb, Bc),
63 StdExpansion3D(Ba.GetNumModes() * Bb.GetNumModes() * Bc.GetNumModes(), Ba,
64 Bb, Bc),
65 StdHexExp(Ba, Bb, Bc), Expansion(geom), Expansion3D(geom),
67 std::bind(&Expansion3D::CreateMatrix, this, std::placeholders::_1),
68 std::string("HexExpMatrix")),
70 this, std::placeholders::_1),
71 std::string("HexExpStaticCondMatrix"))
72{
73}
Expansion3D(SpatialDomains::Geometry3DSharedPtr pGeom)
Definition: Expansion3D.h:59
DNekScalMatSharedPtr CreateMatrix(const MatrixKey &mkey)
DNekScalBlkMatSharedPtr CreateStaticCondMatrix(const MatrixKey &mkey)
Definition: Expansion.cpp:272
Expansion(SpatialDomains::GeometrySharedPtr pGeom)
Definition: Expansion.cpp:43
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: HexExp.h:244
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: HexExp.h:246
StdExpansion()
Default Constructor.
StdHexExp(const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
Definition: StdHexExp.cpp:45

◆ HexExp() [2/2]

Nektar::LocalRegions::HexExp::HexExp ( const HexExp T)

Copy Constructor.

Parameters
THexExp to copy.

Definition at line 80 of file HexExp.cpp.

82 Expansion3D(T), m_matrixManager(T.m_matrixManager),
83 m_staticCondMatrixManager(T.m_staticCondMatrixManager)
84{
85}

◆ ~HexExp()

Nektar::LocalRegions::HexExp::~HexExp ( )
overridedefault

Member Function Documentation

◆ v_AlignVectorToCollapsedDir()

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

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 472 of file HexExp.cpp.

475{
476 ASSERTL1((dir == 0) || (dir == 1) || (dir == 2), "Invalid direction.");
477
478 const int nq0 = m_base[0]->GetNumPoints();
479 const int nq1 = m_base[1]->GetNumPoints();
480 const int nq2 = m_base[2]->GetNumPoints();
481 const int nq = nq0 * nq1 * nq2;
482
483 const Array<TwoD, const NekDouble> &df =
484 m_metricinfo->GetDerivFactors(GetPointsKeys());
485
486 Array<OneD, NekDouble> tmp1(nq); // Quad metric
487
488 Array<OneD, NekDouble> tmp2 = outarray[0]; // Dir1 metric
489 Array<OneD, NekDouble> tmp3 = outarray[1]; // Dir2 metric
490 Array<OneD, NekDouble> tmp4 = outarray[2];
491
492 Vmath::Vcopy(nq, inarray, 1, tmp1, 1); // Dir3 metric
493
494 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
495 {
496 Vmath::Vmul(nq, &df[3 * dir][0], 1, tmp1.get(), 1, tmp2.get(), 1);
497 Vmath::Vmul(nq, &df[3 * dir + 1][0], 1, tmp1.get(), 1, tmp3.get(), 1);
498 Vmath::Vmul(nq, &df[3 * dir + 2][0], 1, tmp1.get(), 1, tmp4.get(), 1);
499 }
500 else
501 {
502 Vmath::Smul(nq, df[3 * dir][0], tmp1.get(), 1, tmp2.get(), 1);
503 Vmath::Smul(nq, df[3 * dir + 1][0], tmp1.get(), 1, tmp3.get(), 1);
504 Vmath::Smul(nq, df[3 * dir + 2][0], tmp1.get(), 1, tmp4.get(), 1);
505 }
506}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:274
const LibUtilities::PointsKeyVector GetPointsKeys() const
Array< OneD, LibUtilities::BasisSharedPtr > m_base
@ eDeformed
Geometry is curved or has non-constant factors.
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.hpp:72
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.hpp:100
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

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

Referenced by v_IProductWRTDerivBase_SumFac().

◆ v_ComputeLaplacianMetric()

void Nektar::LocalRegions::HexExp::v_ComputeLaplacianMetric ( )
overrideprotectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1426 of file HexExp.cpp.

1427{
1428 if (m_metrics.count(eMetricQuadrature) == 0)
1429 {
1431 }
1432
1433 const SpatialDomains::GeomType type = m_metricinfo->GetGtype();
1434 const unsigned int nqtot = GetTotPoints();
1435 const unsigned int dim = 3;
1436 const MetricType m[3][3] = {
1440
1441 for (unsigned int i = 0; i < dim; ++i)
1442 {
1443 for (unsigned int j = i; j < dim; ++j)
1444 {
1445 m_metrics[m[i][j]] = Array<OneD, NekDouble>(nqtot);
1446 const Array<TwoD, const NekDouble> &gmat =
1447 m_metricinfo->GetGmat(GetPointsKeys());
1448 if (type == SpatialDomains::eDeformed)
1449 {
1450 Vmath::Vcopy(nqtot, &gmat[i * dim + j][0], 1,
1451 &m_metrics[m[i][j]][0], 1);
1452 }
1453 else
1454 {
1455 Vmath::Fill(nqtot, gmat[i * dim + j][0], &m_metrics[m[i][j]][0],
1456 1);
1457 }
1459 }
1460 }
1461}
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:134
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:723
GeomType
Indicates the type of element geometry.
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.hpp:54

References Nektar::LocalRegions::Expansion::ComputeQuadratureMetric(), Nektar::SpatialDomains::eDeformed, Nektar::LocalRegions::eMetricLaplacian00, Nektar::LocalRegions::eMetricLaplacian01, Nektar::LocalRegions::eMetricLaplacian02, Nektar::LocalRegions::eMetricLaplacian11, Nektar::LocalRegions::eMetricLaplacian12, Nektar::LocalRegions::eMetricLaplacian22, Nektar::LocalRegions::eMetricQuadrature, Vmath::Fill(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::GetTotPoints(), Nektar::LocalRegions::Expansion::m_metricinfo, Nektar::LocalRegions::Expansion::m_metrics, Nektar::StdRegions::StdExpansion::MultiplyByQuadratureMetric(), and Vmath::Vcopy().

◆ v_ComputeTraceNormal()

void Nektar::LocalRegions::HexExp::v_ComputeTraceNormal ( const int  face)
overrideprotectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 880 of file HexExp.cpp.

881{
882 int i;
883 const SpatialDomains::GeomFactorsSharedPtr &geomFactors =
884 GetGeom()->GetMetricInfo();
885 SpatialDomains::GeomType type = geomFactors->GetGtype();
886
888 for (i = 0; i < ptsKeys.size(); ++i)
889 {
890 // Need at least 2 points for computing normals
891 if (ptsKeys[i].GetNumPoints() == 1)
892 {
893 LibUtilities::PointsKey pKey(2, ptsKeys[i].GetPointsType());
894 ptsKeys[i] = pKey;
895 }
896 }
897
898 const Array<TwoD, const NekDouble> &df =
899 geomFactors->GetDerivFactors(ptsKeys);
900 const Array<OneD, const NekDouble> &jac = geomFactors->GetJac(ptsKeys);
901
902 LibUtilities::BasisKey tobasis0 = GetTraceBasisKey(face, 0);
903 LibUtilities::BasisKey tobasis1 = GetTraceBasisKey(face, 1);
904
905 // Number of quadrature points in face expansion.
906 int nq_face = tobasis0.GetNumPoints() * tobasis1.GetNumPoints();
907
908 int vCoordDim = GetCoordim();
909
910 m_traceNormals[face] = Array<OneD, Array<OneD, NekDouble>>(vCoordDim);
911 Array<OneD, Array<OneD, NekDouble>> &normal = m_traceNormals[face];
912 for (i = 0; i < vCoordDim; ++i)
913 {
914 normal[i] = Array<OneD, NekDouble>(nq_face);
915 }
916
917 size_t nqb = nq_face;
918 size_t nbnd = face;
919 m_elmtBndNormDirElmtLen[nbnd] = Array<OneD, NekDouble>{nqb, 0.0};
920 Array<OneD, NekDouble> &length = m_elmtBndNormDirElmtLen[nbnd];
921
922 // Regular geometry case
923 if ((type == SpatialDomains::eRegular) ||
925 {
926 NekDouble fac;
927 // Set up normals
928 switch (face)
929 {
930 case 0:
931 for (i = 0; i < vCoordDim; ++i)
932 {
933 normal[i][0] = -df[3 * i + 2][0];
934 }
935 break;
936 case 1:
937 for (i = 0; i < vCoordDim; ++i)
938 {
939 normal[i][0] = -df[3 * i + 1][0];
940 }
941 break;
942 case 2:
943 for (i = 0; i < vCoordDim; ++i)
944 {
945 normal[i][0] = df[3 * i][0];
946 }
947 break;
948 case 3:
949 for (i = 0; i < vCoordDim; ++i)
950 {
951 normal[i][0] = df[3 * i + 1][0];
952 }
953 break;
954 case 4:
955 for (i = 0; i < vCoordDim; ++i)
956 {
957 normal[i][0] = -df[3 * i][0];
958 }
959 break;
960 case 5:
961 for (i = 0; i < vCoordDim; ++i)
962 {
963 normal[i][0] = df[3 * i + 2][0];
964 }
965 break;
966 default:
967 ASSERTL0(false, "face is out of range (edge < 5)");
968 }
969
970 // normalise
971 fac = 0.0;
972 for (i = 0; i < vCoordDim; ++i)
973 {
974 fac += normal[i][0] * normal[i][0];
975 }
976 fac = 1.0 / sqrt(fac);
977
978 Vmath::Fill(nqb, fac, length, 1);
979 for (i = 0; i < vCoordDim; ++i)
980 {
981 Vmath::Fill(nq_face, fac * normal[i][0], normal[i], 1);
982 }
983 }
984 else // Set up deformed normals
985 {
986 int j, k;
987
988 int nqe0 = ptsKeys[0].GetNumPoints();
989 int nqe1 = ptsKeys[1].GetNumPoints();
990 int nqe2 = ptsKeys[2].GetNumPoints();
991 int nqe01 = nqe0 * nqe1;
992 int nqe02 = nqe0 * nqe2;
993 int nqe12 = nqe1 * nqe2;
994
995 int nqe;
996 if (face == 0 || face == 5)
997 {
998 nqe = nqe01;
999 }
1000 else if (face == 1 || face == 3)
1001 {
1002 nqe = nqe02;
1003 }
1004 else
1005 {
1006 nqe = nqe12;
1007 }
1008
1009 LibUtilities::PointsKey points0;
1010 LibUtilities::PointsKey points1;
1011
1012 Array<OneD, NekDouble> faceJac(nqe);
1013 Array<OneD, NekDouble> normals(vCoordDim * nqe, 0.0);
1014
1015 // Extract Jacobian along face and recover local
1016 // derivates (dx/dr) for polynomial interpolation by
1017 // multiplying m_gmat by jacobian
1018 switch (face)
1019 {
1020 case 0:
1021 for (j = 0; j < nqe; ++j)
1022 {
1023 normals[j] = -df[2][j] * jac[j];
1024 normals[nqe + j] = -df[5][j] * jac[j];
1025 normals[2 * nqe + j] = -df[8][j] * jac[j];
1026 faceJac[j] = jac[j];
1027 }
1028
1029 points0 = ptsKeys[0];
1030 points1 = ptsKeys[1];
1031 break;
1032 case 1:
1033 for (j = 0; j < nqe0; ++j)
1034 {
1035 for (k = 0; k < nqe2; ++k)
1036 {
1037 int idx = j + nqe01 * k;
1038 normals[j + k * nqe0] = -df[1][idx] * jac[idx];
1039 normals[nqe + j + k * nqe0] = -df[4][idx] * jac[idx];
1040 normals[2 * nqe + j + k * nqe0] =
1041 -df[7][idx] * jac[idx];
1042 faceJac[j + k * nqe0] = jac[idx];
1043 }
1044 }
1045 points0 = ptsKeys[0];
1046 points1 = ptsKeys[2];
1047 break;
1048 case 2:
1049 for (j = 0; j < nqe1; ++j)
1050 {
1051 for (k = 0; k < nqe2; ++k)
1052 {
1053 int idx = nqe0 - 1 + nqe0 * j + nqe01 * k;
1054 normals[j + k * nqe1] = df[0][idx] * jac[idx];
1055 normals[nqe + j + k * nqe1] = df[3][idx] * jac[idx];
1056 normals[2 * nqe + j + k * nqe1] = df[6][idx] * jac[idx];
1057 faceJac[j + k * nqe1] = jac[idx];
1058 }
1059 }
1060 points0 = ptsKeys[1];
1061 points1 = ptsKeys[2];
1062 break;
1063 case 3:
1064 for (j = 0; j < nqe0; ++j)
1065 {
1066 for (k = 0; k < nqe2; ++k)
1067 {
1068 int idx = nqe0 * (nqe1 - 1) + j + nqe01 * k;
1069 normals[j + k * nqe0] = df[1][idx] * jac[idx];
1070 normals[nqe + j + k * nqe0] = df[4][idx] * jac[idx];
1071 normals[2 * nqe + j + k * nqe0] = df[7][idx] * jac[idx];
1072 faceJac[j + k * nqe0] = jac[idx];
1073 }
1074 }
1075 points0 = ptsKeys[0];
1076 points1 = ptsKeys[2];
1077 break;
1078 case 4:
1079 for (j = 0; j < nqe1; ++j)
1080 {
1081 for (k = 0; k < nqe2; ++k)
1082 {
1083 int idx = j * nqe0 + nqe01 * k;
1084 normals[j + k * nqe1] = -df[0][idx] * jac[idx];
1085 normals[nqe + j + k * nqe1] = -df[3][idx] * jac[idx];
1086 normals[2 * nqe + j + k * nqe1] =
1087 -df[6][idx] * jac[idx];
1088 faceJac[j + k * nqe1] = jac[idx];
1089 }
1090 }
1091 points0 = ptsKeys[1];
1092 points1 = ptsKeys[2];
1093 break;
1094 case 5:
1095 for (j = 0; j < nqe01; ++j)
1096 {
1097 int idx = j + nqe01 * (nqe2 - 1);
1098 normals[j] = df[2][idx] * jac[idx];
1099 normals[nqe + j] = df[5][idx] * jac[idx];
1100 normals[2 * nqe + j] = df[8][idx] * jac[idx];
1101 faceJac[j] = jac[idx];
1102 }
1103 points0 = ptsKeys[0];
1104 points1 = ptsKeys[1];
1105 break;
1106 default:
1107 ASSERTL0(false, "face is out of range (face < 5)");
1108 }
1109
1110 Array<OneD, NekDouble> work(nq_face, 0.0);
1111 // Interpolate Jacobian and invert
1112 LibUtilities::Interp2D(points0, points1, faceJac,
1113 tobasis0.GetPointsKey(), tobasis1.GetPointsKey(),
1114 work);
1115
1116 Vmath::Sdiv(nq_face, 1.0, &work[0], 1, &work[0], 1);
1117
1118 // interpolate
1119 for (i = 0; i < GetCoordim(); ++i)
1120 {
1121 LibUtilities::Interp2D(points0, points1, &normals[i * nqe],
1122 tobasis0.GetPointsKey(),
1123 tobasis1.GetPointsKey(), &normal[i][0]);
1124 Vmath::Vmul(nq_face, work, 1, normal[i], 1, normal[i], 1);
1125 }
1126
1127 // normalise normal vectors
1128 Vmath::Zero(nq_face, work, 1);
1129 for (i = 0; i < GetCoordim(); ++i)
1130 {
1131 Vmath::Vvtvp(nq_face, normal[i], 1, normal[i], 1, work, 1, work, 1);
1132 }
1133
1134 Vmath::Vsqrt(nq_face, work, 1, work, 1);
1135 Vmath::Sdiv(nq_face, 1.0, work, 1, work, 1);
1136
1137 Vmath::Vcopy(nqb, work, 1, length, 1);
1138
1139 for (i = 0; i < GetCoordim(); ++i)
1140 {
1141 Vmath::Vmul(nq_face, normal[i], 1, work, 1, normal[i], 1);
1142 }
1143 }
1144}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
std::map< int, NormalVector > m_traceNormals
Definition: Expansion.h:276
std::map< int, Array< OneD, NekDouble > > m_elmtBndNormDirElmtLen
the element length in each element boundary(Vertex, edge or face) normal direction calculated based o...
Definition: Expansion.h:286
SpatialDomains::GeometrySharedPtr GetGeom() const
Definition: Expansion.cpp:167
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:299
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:205
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:218
void Interp2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
this function interpolates a 2D function evaluated at the quadrature points of the 2D basis,...
Definition: Interp.cpp:101
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:231
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition: GeomFactors.h:60
@ eRegular
Geometry is straight-sided with constant geometric factors.
@ eMovingRegular
Currently unused.
double NekDouble
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
Definition: Vmath.hpp:340
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.hpp:366
void Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/x.
Definition: Vmath.hpp:154
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

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::HexExp::v_CreateStdMatrix ( const StdRegions::StdMatrixKey mkey)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdHexExp.

Definition at line 1321 of file HexExp.cpp.

1322{
1323 LibUtilities::BasisKey bkey0 = m_base[0]->GetBasisKey();
1324 LibUtilities::BasisKey bkey1 = m_base[1]->GetBasisKey();
1325 LibUtilities::BasisKey bkey2 = m_base[2]->GetBasisKey();
1326
1329
1330 return tmp->GetStdMatrix(mkey);
1331}
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< StdHexExp > StdHexExpSharedPtr
Definition: StdHexExp.h:228

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

◆ v_DetShapeType()

LibUtilities::ShapeType Nektar::LocalRegions::HexExp::v_DetShapeType ( ) const
overrideprotectedvirtual

Return the region shape using the enum-list of ShapeType.

Reimplemented from Nektar::StdRegions::StdHexExp.

Definition at line 657 of file HexExp.cpp.

658{
660}

References Nektar::LibUtilities::eHexahedron.

◆ v_DropLocMatrix()

void Nektar::LocalRegions::HexExp::v_DropLocMatrix ( const MatrixKey mkey)
overrideprotectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1338 of file HexExp.cpp.

1339{
1340 m_matrixManager.DeleteObject(mkey);
1341}

References m_matrixManager.

◆ v_DropLocStaticCondMatrix()

void Nektar::LocalRegions::HexExp::v_DropLocStaticCondMatrix ( const MatrixKey mkey)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1348 of file HexExp.cpp.

1349{
1350 m_staticCondMatrixManager.DeleteObject(mkey);
1351}

References m_staticCondMatrixManager.

◆ v_ExtractDataToCoeffs()

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

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 662 of file HexExp.cpp.

666{
667 int data_order0 = nummodes[mode_offset];
668 int fillorder0 = min(m_base[0]->GetNumModes(), data_order0);
669 int data_order1 = nummodes[mode_offset + 1];
670 int order1 = m_base[1]->GetNumModes();
671 int fillorder1 = min(order1, data_order1);
672 int data_order2 = nummodes[mode_offset + 2];
673 int order2 = m_base[2]->GetNumModes();
674 int fillorder2 = min(order2, data_order2);
675
676 // Check if same basis
677 if (fromType[0] != m_base[0]->GetBasisType() ||
678 fromType[1] != m_base[1]->GetBasisType() ||
679 fromType[2] != m_base[2]->GetBasisType())
680 {
681 // Construct a hex with the appropriate basis type at our
682 // quadrature points, and one more to do a forwards
683 // transform. We can then copy the output to coeffs.
684 StdRegions::StdHexExp tmpHex(
685 LibUtilities::BasisKey(fromType[0], data_order0,
686 m_base[0]->GetPointsKey()),
687 LibUtilities::BasisKey(fromType[1], data_order1,
688 m_base[1]->GetPointsKey()),
689 LibUtilities::BasisKey(fromType[2], data_order2,
690 m_base[2]->GetPointsKey()));
691 StdRegions::StdHexExp tmpHex2(m_base[0]->GetBasisKey(),
692 m_base[1]->GetBasisKey(),
693 m_base[2]->GetBasisKey());
694
695 Array<OneD, const NekDouble> tmpData(tmpHex.GetNcoeffs(), data);
696 Array<OneD, NekDouble> tmpBwd(tmpHex2.GetTotPoints());
697 Array<OneD, NekDouble> tmpOut(tmpHex2.GetNcoeffs());
698
699 tmpHex.BwdTrans(tmpData, tmpBwd);
700 tmpHex2.FwdTrans(tmpBwd, tmpOut);
701 Vmath::Vcopy(tmpOut.size(), &tmpOut[0], 1, coeffs, 1);
702
703 return;
704 }
705
706 switch (m_base[0]->GetBasisType())
707 {
709 {
710 int i, j;
711 int cnt = 0;
712 int cnt1 = 0;
713
715 "Extraction routine not set up for this basis");
717 "Extraction routine not set up for this basis");
718
719 Vmath::Zero(m_ncoeffs, coeffs, 1);
720 for (j = 0; j < fillorder0; ++j)
721 {
722 for (i = 0; i < fillorder1; ++i)
723 {
724 Vmath::Vcopy(fillorder2, &data[cnt], 1, &coeffs[cnt1], 1);
725 cnt += data_order2;
726 cnt1 += order2;
727 }
728
729 // count out data for j iteration
730 for (i = fillorder1; i < data_order1; ++i)
731 {
732 cnt += data_order2;
733 }
734
735 for (i = fillorder1; i < order1; ++i)
736 {
737 cnt1 += order2;
738 }
739 }
740 break;
741 }
743 {
744 LibUtilities::PointsKey p0(nummodes[0],
746 LibUtilities::PointsKey p1(nummodes[1],
748 LibUtilities::PointsKey p2(nummodes[2],
750 LibUtilities::PointsKey t0(m_base[0]->GetNumModes(),
752 LibUtilities::PointsKey t1(m_base[1]->GetNumModes(),
754 LibUtilities::PointsKey t2(m_base[2]->GetNumModes(),
756 LibUtilities::Interp3D(p0, p1, p2, data, t0, t1, t2, coeffs);
757 }
758 break;
759 default:
760 ASSERTL0(false, "basis is either not set up or not "
761 "hierarchicial");
762 }
763}
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:156
void Interp3D(const BasisKey &fbasis0, const BasisKey &fbasis1, const BasisKey &fbasis2, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, const BasisKey &tbasis2, Array< OneD, NekDouble > &to)
this function interpolates a 3D function evaluated at the quadrature points of the 3D basis,...
Definition: Interp.cpp:162
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:56
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48

References ASSERTL0, ASSERTL1, Nektar::StdRegions::StdExpansion::BwdTrans(), Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::FwdTrans(), Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), Nektar::StdRegions::StdExpansion::GetTotPoints(), Nektar::LibUtilities::Interp3D(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Vmath::Vcopy(), and Vmath::Zero().

◆ v_FwdTrans()

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

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

Parameters
inarrayInput array
outarrayOutput array

Reimplemented from Nektar::StdRegions::StdHexExp.

Definition at line 297 of file HexExp.cpp.

299{
300 if (m_base[0]->Collocation() && m_base[1]->Collocation() &&
301 m_base[2]->Collocation())
302 {
303 Vmath::Vcopy(GetNcoeffs(), &inarray[0], 1, &outarray[0], 1);
304 }
305 else
306 {
307 IProductWRTBase(inarray, outarray);
308
309 // get Mass matrix inverse
310 MatrixKey masskey(StdRegions::eInvMass, DetShapeType(), *this);
311 DNekScalMatSharedPtr matsys = m_matrixManager[masskey];
312
313 // copy inarray in case inarray == outarray
314 DNekVec in(m_ncoeffs, outarray);
315 DNekVec out(m_ncoeffs, outarray, eWrapper);
316
317 out = (*matsys) * in;
318 }
319}
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:124
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:528
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:367
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
NekVector< NekDouble > DNekVec
Definition: NekTypeDefs.hpp:48

References Nektar::StdRegions::StdExpansion::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::HexExp::v_GenMatrix ( const StdRegions::StdMatrixKey mkey)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdHexExp.

Definition at line 1299 of file HexExp.cpp.

1300{
1301 DNekMatSharedPtr returnval;
1302
1303 switch (mkey.GetMatrixType())
1304 {
1312 returnval = Expansion3D::v_GenMatrix(mkey);
1313 break;
1314 default:
1315 returnval = StdHexExp::v_GenMatrix(mkey);
1316 }
1317
1318 return returnval;
1319}
DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey) override
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::HexExp::v_GetCoord ( const Array< OneD, const NekDouble > &  Lcoords,
Array< OneD, NekDouble > &  coords 
)
overrideprotectedvirtual

Retrieves the physical coordinates of a given set of reference coordinates.

Parameters
LcoordsLocal coordinates in reference space.
coordsCorresponding coordinates in physical space.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 628 of file HexExp.cpp.

630{
631 int i;
632
633 ASSERTL1(Lcoords[0] >= -1.0 && Lcoords[0] <= 1.0 && Lcoords[1] >= -1.0 &&
634 Lcoords[1] <= 1.0 && Lcoords[2] >= -1.0 && Lcoords[2] <= 1.0,
635 "Local coordinates are not in region [-1,1]");
636
637 m_geom->FillGeom();
638
639 for (i = 0; i < m_geom->GetCoordim(); ++i)
640 {
641 coords[i] = m_geom->GetCoord(i, Lcoords);
642 }
643}
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:273

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

◆ v_GetCoords()

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

Reimplemented from Nektar::StdRegions::StdHexExp.

Definition at line 645 of file HexExp.cpp.

648{
649 Expansion::v_GetCoords(coords_0, coords_1, coords_2);
650}
void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3) override
Definition: Expansion.cpp:530

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

◆ v_GetLinStdExp()

StdRegions::StdExpansionSharedPtr Nektar::LocalRegions::HexExp::v_GetLinStdExp ( void  ) const
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 608 of file HexExp.cpp.

609{
610 LibUtilities::BasisKey bkey0(m_base[0]->GetBasisType(), 2,
611 m_base[0]->GetPointsKey());
612 LibUtilities::BasisKey bkey1(m_base[1]->GetBasisType(), 2,
613 m_base[1]->GetPointsKey());
614 LibUtilities::BasisKey bkey2(m_base[2]->GetBasisType(), 2,
615 m_base[2]->GetPointsKey());
616
618 bkey2);
619}

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

◆ v_GetLocMatrix()

DNekScalMatSharedPtr Nektar::LocalRegions::HexExp::v_GetLocMatrix ( const MatrixKey mkey)
overrideprotectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1333 of file HexExp.cpp.

1334{
1335 return m_matrixManager[mkey];
1336}

References m_matrixManager.

◆ v_GetLocStaticCondMatrix()

DNekScalBlkMatSharedPtr Nektar::LocalRegions::HexExp::v_GetLocStaticCondMatrix ( const MatrixKey mkey)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1343 of file HexExp.cpp.

1344{
1345 return m_staticCondMatrixManager[mkey];
1346}

References m_staticCondMatrixManager.

◆ v_GetStdExp()

StdRegions::StdExpansionSharedPtr Nektar::LocalRegions::HexExp::v_GetStdExp ( void  ) const
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 601 of file HexExp.cpp.

602{
604 m_base[0]->GetBasisKey(), m_base[1]->GetBasisKey(),
605 m_base[2]->GetBasisKey());
606}

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

◆ v_GetTracePhysMap()

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

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 765 of file HexExp.cpp.

766{
767 int nquad0 = m_base[0]->GetNumPoints();
768 int nquad1 = m_base[1]->GetNumPoints();
769 int nquad2 = m_base[2]->GetNumPoints();
770
771 int nq0 = 0;
772 int nq1 = 0;
773
774 switch (face)
775 {
776 case 0:
777 nq0 = nquad0;
778 nq1 = nquad1;
779
780 // Directions A and B positive
781 if (outarray.size() != nq0 * nq1)
782 {
783 outarray = Array<OneD, int>(nq0 * nq1);
784 }
785
786 for (int i = 0; i < nquad0 * nquad1; ++i)
787 {
788 outarray[i] = i;
789 }
790
791 break;
792 case 1:
793 nq0 = nquad0;
794 nq1 = nquad2;
795 // Direction A and B positive
796 if (outarray.size() != nq0 * nq1)
797 {
798 outarray = Array<OneD, int>(nq0 * nq1);
799 }
800
801 // Direction A and B positive
802 for (int k = 0; k < nquad2; k++)
803 {
804 for (int i = 0; i < nquad0; ++i)
805 {
806 outarray[k * nquad0 + i] = nquad0 * nquad1 * k + i;
807 }
808 }
809 break;
810 case 2:
811 nq0 = nquad1;
812 nq1 = nquad2;
813
814 // Direction A and B positive
815 if (outarray.size() != nq0 * nq1)
816 {
817 outarray = Array<OneD, int>(nq0 * nq1);
818 }
819
820 for (int i = 0; i < nquad1 * nquad2; i++)
821 {
822 outarray[i] = nquad0 - 1 + i * nquad0;
823 }
824 break;
825 case 3:
826 nq0 = nquad0;
827 nq1 = nquad2;
828
829 // Direction A and B positive
830 if (outarray.size() != nq0 * nq1)
831 {
832 outarray = Array<OneD, int>(nq0 * nq1);
833 }
834
835 for (int k = 0; k < nquad2; k++)
836 {
837 for (int i = 0; i < nquad0; i++)
838 {
839 outarray[k * nquad0 + i] =
840 (nquad0 * (nquad1 - 1)) + (k * nquad0 * nquad1) + i;
841 }
842 }
843 break;
844 case 4:
845 nq0 = nquad1;
846 nq1 = nquad2;
847
848 // Direction A and B positive
849 if (outarray.size() != nq0 * nq1)
850 {
851 outarray = Array<OneD, int>(nq0 * nq1);
852 }
853
854 for (int i = 0; i < nquad1 * nquad2; i++)
855 {
856 outarray[i] = i * nquad0;
857 }
858 break;
859 case 5:
860 nq0 = nquad0;
861 nq1 = nquad1;
862 // Directions A and B positive
863 if (outarray.size() != nq0 * nq1)
864 {
865 outarray = Array<OneD, int>(nq0 * nq1);
866 }
867
868 for (int i = 0; i < nquad0 * nquad1; i++)
869 {
870 outarray[i] = nquad0 * nquad1 * (nquad2 - 1) + i;
871 }
872
873 break;
874 default:
875 ASSERTL0(false, "face value (> 5) is out of range");
876 break;
877 }
878}

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

◆ v_HelmholtzMatrixOp()

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

Reimplemented from Nektar::StdRegions::StdHexExp.

Definition at line 1193 of file HexExp.cpp.

1196{
1197 HexExp::v_HelmholtzMatrixOp_MatFree(inarray, outarray, mkey);
1198}
void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override

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

◆ v_Integral()

NekDouble Nektar::LocalRegions::HexExp::v_Integral ( const Array< OneD, const NekDouble > &  inarray)
overrideprotectedvirtual

Integrate the physical point list inarray over region.

Parameters
inarraydefinition of function to be returned at quadrature points 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 101 of file HexExp.cpp.

102{
103 int nquad0 = m_base[0]->GetNumPoints();
104 int nquad1 = m_base[1]->GetNumPoints();
105 int nquad2 = m_base[2]->GetNumPoints();
106 Array<OneD, const NekDouble> jac = m_metricinfo->GetJac(GetPointsKeys());
107 NekDouble returnVal;
108 Array<OneD, NekDouble> tmp(nquad0 * nquad1 * nquad2);
109
110 // multiply inarray with Jacobian
111
112 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
113 {
114 Vmath::Vmul(nquad0 * nquad1 * nquad2, &jac[0], 1,
115 (NekDouble *)&inarray[0], 1, &tmp[0], 1);
116 }
117 else
118 {
119 Vmath::Smul(nquad0 * nquad1 * nquad2, (NekDouble)jac[0],
120 (NekDouble *)&inarray[0], 1, &tmp[0], 1);
121 }
122
123 // call StdHexExp version;
124 returnVal = StdHexExp::v_Integral(tmp);
125
126 return returnVal;
127}

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::HexExp::v_IProductWRTBase ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Calculate the inner product of inarray with respect to the elements basis.

Parameters
inarrayInput array of physical space data.
outarrayOutput array of data.

Reimplemented from Nektar::StdRegions::StdHexExp.

Definition at line 332 of file HexExp.cpp.

334{
335 HexExp::v_IProductWRTBase_SumFac(inarray, outarray);
336}
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
Calculate the inner product of inarray with respect to the given basis B = base0 * base1 * base2.
Definition: HexExp.cpp:370

References v_IProductWRTBase_SumFac().

◆ v_IProductWRTBase_SumFac()

void Nektar::LocalRegions::HexExp::v_IProductWRTBase_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
bool  multiplybyweights = true 
)
overrideprotectedvirtual

Calculate the inner product of inarray with respect to the given basis B = base0 * base1 * base2.

\( \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} (\xi_{1i}) \psi_{q}^{a} (\xi_{2j}) \psi_{r}^{a} (\xi_{3k}) w_i w_j w_k u(\xi_{1,i} \xi_{2,j} \xi_{3,k}) J_{i,j,k}\\ & = & \sum_{i=0}^{nq_0} \psi_p^a(\xi_{1,i}) \sum_{j=0}^{nq_1} \psi_{q}^a(\xi_{2,j}) \sum_{k=0}^{nq_2} \psi_{r}^a u(\xi_{1i},\xi_{2j},\xi_{3k}) J_{i,j,k} \end{array} \)
where \( \phi_{pqr} (\xi_1 , \xi_2 , \xi_3) = \psi_p^a ( \xi_1) \psi_{q}^a (\xi_2) \psi_{r}^a (\xi_3) \)
which can be implemented as
\(f_{r} (\xi_{3k}) = \sum_{k=0}^{nq_3} \psi_{r}^a u(\xi_{1i},\xi_{2j},\xi_{3k}) J_{i,j,k} = {\bf B_3 U} \)
\( g_{q} (\xi_{3k}) = \sum_{j=0}^{nq_1} \psi_{q}^a (\xi_{2j}) f_{r} (\xi_{3k}) = {\bf B_2 F} \)
\( (\phi_{pqr}, u)_{\delta} = \sum_{k=0}^{nq_0} \psi_{p}^a (\xi_{3k}) g_{q} (\xi_{3k}) = {\bf B_1 G} \)

Parameters
base0Basis to integrate wrt in first dimension.
base1Basis to integrate wrt in second dimension.
base2Basis to integrate wrt in third dimension.
inarrayInput array.
outarrayOutput array.
coll_check(not used)

Reimplemented from Nektar::StdRegions::StdHexExp.

Definition at line 370 of file HexExp.cpp.

373{
374 int nquad0 = m_base[0]->GetNumPoints();
375 int nquad1 = m_base[1]->GetNumPoints();
376 int nquad2 = m_base[2]->GetNumPoints();
377 int order0 = m_base[0]->GetNumModes();
378 int order1 = m_base[1]->GetNumModes();
379
380 Array<OneD, NekDouble> wsp(nquad0 * nquad1 * (nquad2 + order0) +
381 order0 * order1 * nquad2);
382
383 if (multiplybyweights)
384 {
385 Array<OneD, NekDouble> tmp(inarray.size());
386
387 MultiplyByQuadratureMetric(inarray, tmp);
389 m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
390 tmp, outarray, wsp, true, true, true);
391 }
392 else
393 {
395 m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
396 inarray, outarray, wsp, true, true, true);
397 }
398}
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::HexExp::v_IProductWRTDerivBase ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdHexExp.

Definition at line 400 of file HexExp.cpp.

403{
404 HexExp::v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
405}
void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Calculates the inner product .
Definition: HexExp.cpp:427

References v_IProductWRTDerivBase_SumFac().

◆ v_IProductWRTDerivBase_SumFac()

void Nektar::LocalRegions::HexExp::v_IProductWRTDerivBase_SumFac ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

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}(\xi_1,\xi_2,\xi_3) = \phi_1(\xi_1)\phi_2(\xi_2)\phi_3(\xi_3)\), in the hexahedral element, this is straightforward and yields the result

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

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

Reimplemented from Nektar::StdRegions::StdHexExp.

Definition at line 427 of file HexExp.cpp.

430{
431 ASSERTL1((dir == 0) || (dir == 1) || (dir == 2), "Invalid direction.");
432
433 const int nq0 = m_base[0]->GetNumPoints();
434 const int nq1 = m_base[1]->GetNumPoints();
435 const int nq2 = m_base[2]->GetNumPoints();
436 const int nq = nq0 * nq1 * nq2;
437 const int nm0 = m_base[0]->GetNumModes();
438 const int nm1 = m_base[1]->GetNumModes();
439
440 Array<OneD, NekDouble> alloc(4 * nq + m_ncoeffs + nm0 * nq2 * (nq1 + nm1));
441 Array<OneD, NekDouble> tmp1(alloc); // Quad metric
442 Array<OneD, NekDouble> tmp2(alloc + nq); // Dir1 metric
443 Array<OneD, NekDouble> tmp3(alloc + 2 * nq); // Dir2 metric
444 Array<OneD, NekDouble> tmp4(alloc + 3 * nq); // Dir3 metric
445 Array<OneD, NekDouble> tmp5(alloc + 4 * nq); // iprod tmp
446 Array<OneD, NekDouble> wsp(tmp5 + m_ncoeffs); // Wsp
447
448 MultiplyByQuadratureMetric(inarray, tmp1);
449
450 Array<OneD, Array<OneD, NekDouble>> tmp2D{3};
451 tmp2D[0] = tmp2;
452 tmp2D[1] = tmp3;
453 tmp2D[2] = tmp4;
454
455 HexExp::v_AlignVectorToCollapsedDir(dir, tmp1, tmp2D);
456
457 IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
458 m_base[2]->GetBdata(), tmp2, outarray, wsp,
459 false, true, true);
460
461 IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
462 m_base[2]->GetBdata(), tmp3, tmp5, wsp, true,
463 false, true);
464 Vmath::Vadd(m_ncoeffs, tmp5, 1, outarray, 1, outarray, 1);
465
466 IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetBdata(),
467 m_base[2]->GetDbdata(), tmp4, tmp5, wsp, true,
468 true, false);
469 Vmath::Vadd(m_ncoeffs, tmp5, 1, outarray, 1, outarray, 1);
470}
void v_AlignVectorToCollapsedDir(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
Definition: HexExp.cpp:472
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.hpp:180

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

Referenced by v_IProductWRTDerivBase().

◆ v_IProductWRTDirectionalDerivBase()

void Nektar::LocalRegions::HexExp::v_IProductWRTDirectionalDerivBase ( const Array< OneD, const NekDouble > &  direction,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
inlineoverrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 115 of file HexExp.h.

119 {
120 IProductWRTDirectionalDerivBase_SumFac(direction, inarray, outarray);
121 }
void IProductWRTDirectionalDerivBase_SumFac(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)

References Nektar::StdRegions::StdExpansion::IProductWRTDirectionalDerivBase_SumFac().

◆ v_IProductWRTDirectionalDerivBase_SumFac()

void Nektar::LocalRegions::HexExp::v_IProductWRTDirectionalDerivBase_SumFac ( const Array< OneD, const NekDouble > &  direction,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual
Parameters
dirVector direction in which to take the derivative.
inarrayThe function \( u \).
outarrayValue of the inner product.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 514 of file HexExp.cpp.

518{
519 int shapedim = 3;
520 const int nq0 = m_base[0]->GetNumPoints();
521 const int nq1 = m_base[1]->GetNumPoints();
522 const int nq2 = m_base[2]->GetNumPoints();
523 const int nq = nq0 * nq1 * nq2;
524 const int nm0 = m_base[0]->GetNumModes();
525 const int nm1 = m_base[1]->GetNumModes();
526
527 const Array<TwoD, const NekDouble> &df =
528 m_metricinfo->GetDerivFactors(GetPointsKeys());
529
530 Array<OneD, NekDouble> alloc(4 * nq + m_ncoeffs + nm0 * nq2 * (nq1 + nm1));
531 Array<OneD, NekDouble> tmp1(alloc); // Quad metric
532 Array<OneD, NekDouble> tmp2(alloc + nq); // Dir1 metric
533 Array<OneD, NekDouble> tmp3(alloc + 2 * nq); // Dir2 metric
534 Array<OneD, NekDouble> tmp4(alloc + 3 * nq); // Dir3 metric
535 Array<OneD, NekDouble> tmp5(alloc + 4 * nq); // iprod tmp
536 Array<OneD, NekDouble> wsp(tmp5 + m_ncoeffs); // Wsp
537
538 MultiplyByQuadratureMetric(inarray, tmp1);
539
540 Array<OneD, Array<OneD, NekDouble>> dfdir(shapedim);
541 Expansion::ComputeGmatcdotMF(df, direction, dfdir);
542
543 Vmath::Vmul(nq, &dfdir[0][0], 1, tmp1.get(), 1, tmp2.get(), 1);
544 Vmath::Vmul(nq, &dfdir[1][0], 1, tmp1.get(), 1, tmp3.get(), 1);
545 Vmath::Vmul(nq, &dfdir[2][0], 1, tmp1.get(), 1, tmp4.get(), 1);
546
547 IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
548 m_base[2]->GetBdata(), tmp2, outarray, wsp,
549 false, true, true);
550
551 IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
552 m_base[2]->GetBdata(), tmp3, tmp5, wsp, true,
553 false, true);
554
555 Vmath::Vadd(m_ncoeffs, tmp5, 1, outarray, 1, outarray, 1);
556
557 IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetBdata(),
558 m_base[2]->GetDbdata(), tmp4, tmp5, wsp, true,
559 true, false);
560
561 Vmath::Vadd(m_ncoeffs, tmp5, 1, outarray, 1, outarray, 1);
562}
void ComputeGmatcdotMF(const Array< TwoD, const NekDouble > &df, const Array< OneD, const NekDouble > &direction, Array< OneD, Array< OneD, NekDouble > > &dfdir)
Definition: Expansion.cpp:603

References Nektar::LocalRegions::Expansion::ComputeGmatcdotMF(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion3D::IProductWRTBase_SumFacKernel(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_metricinfo, Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::StdRegions::StdExpansion::MultiplyByQuadratureMetric(), Vmath::Vadd(), and Vmath::Vmul().

◆ v_LaplacianMatrixOp() [1/2]

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

Reimplemented from Nektar::StdRegions::StdHexExp.

Definition at line 1156 of file HexExp.cpp.

1159{
1160 HexExp::v_LaplacianMatrixOp_MatFree(inarray, outarray, mkey);
1161}
void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override

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

◆ v_LaplacianMatrixOp() [2/2]

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

Reimplemented from Nektar::StdRegions::StdHexExp.

Definition at line 1163 of file HexExp.cpp.

1167{
1168 StdExpansion::LaplacianMatrixOp_MatFree(k1, k2, inarray, outarray, mkey);
1169}

◆ v_LaplacianMatrixOp_MatFree_Kernel()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1353 of file HexExp.cpp.

1356{
1357 // This implementation is only valid when there are no
1358 // coefficients associated to the Laplacian operator
1359 if (m_metrics.count(eMetricLaplacian00) == 0)
1360 {
1362 }
1363
1364 int nquad0 = m_base[0]->GetNumPoints();
1365 int nquad1 = m_base[1]->GetNumPoints();
1366 int nquad2 = m_base[2]->GetNumPoints();
1367 int nqtot = nquad0 * nquad1 * nquad2;
1368
1369 ASSERTL1(wsp.size() >= 6 * nqtot, "Insufficient workspace size.");
1370
1371 const Array<OneD, const NekDouble> &base0 = m_base[0]->GetBdata();
1372 const Array<OneD, const NekDouble> &base1 = m_base[1]->GetBdata();
1373 const Array<OneD, const NekDouble> &base2 = m_base[2]->GetBdata();
1374 const Array<OneD, const NekDouble> &dbase0 = m_base[0]->GetDbdata();
1375 const Array<OneD, const NekDouble> &dbase1 = m_base[1]->GetDbdata();
1376 const Array<OneD, const NekDouble> &dbase2 = m_base[2]->GetDbdata();
1377 const Array<OneD, const NekDouble> &metric00 =
1379 const Array<OneD, const NekDouble> &metric01 =
1381 const Array<OneD, const NekDouble> &metric02 =
1383 const Array<OneD, const NekDouble> &metric11 =
1385 const Array<OneD, const NekDouble> &metric12 =
1387 const Array<OneD, const NekDouble> &metric22 =
1389
1390 // Allocate temporary storage
1391 Array<OneD, NekDouble> wsp0(wsp);
1392 Array<OneD, NekDouble> wsp1(wsp + 1 * nqtot);
1393 Array<OneD, NekDouble> wsp2(wsp + 2 * nqtot);
1394 Array<OneD, NekDouble> wsp3(wsp + 3 * nqtot);
1395 Array<OneD, NekDouble> wsp4(wsp + 4 * nqtot);
1396 Array<OneD, NekDouble> wsp5(wsp + 5 * nqtot);
1397
1398 StdExpansion3D::PhysTensorDeriv(inarray, wsp0, wsp1, wsp2);
1399
1400 // wsp0 = k = g0 * wsp1 + g1 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
1401 // wsp2 = l = g1 * wsp1 + g2 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
1402 // where g0, g1 and g2 are the metric terms set up in the GeomFactors class
1403 // especially for this purpose
1404 Vmath::Vvtvvtp(nqtot, &metric00[0], 1, &wsp0[0], 1, &metric01[0], 1,
1405 &wsp1[0], 1, &wsp3[0], 1);
1406 Vmath::Vvtvp(nqtot, &metric02[0], 1, &wsp2[0], 1, &wsp3[0], 1, &wsp3[0], 1);
1407 Vmath::Vvtvvtp(nqtot, &metric01[0], 1, &wsp0[0], 1, &metric11[0], 1,
1408 &wsp1[0], 1, &wsp4[0], 1);
1409 Vmath::Vvtvp(nqtot, &metric12[0], 1, &wsp2[0], 1, &wsp4[0], 1, &wsp4[0], 1);
1410 Vmath::Vvtvvtp(nqtot, &metric02[0], 1, &wsp0[0], 1, &metric12[0], 1,
1411 &wsp1[0], 1, &wsp5[0], 1);
1412 Vmath::Vvtvp(nqtot, &metric22[0], 1, &wsp2[0], 1, &wsp5[0], 1, &wsp5[0], 1);
1413
1414 // outarray = m = (D_xi1 * B)^T * k
1415 // wsp1 = n = (D_xi2 * B)^T * l
1416 IProductWRTBase_SumFacKernel(dbase0, base1, base2, wsp3, outarray, wsp0,
1417 false, true, true);
1418 IProductWRTBase_SumFacKernel(base0, dbase1, base2, wsp4, wsp2, wsp0, true,
1419 false, true);
1420 Vmath::Vadd(m_ncoeffs, wsp2.get(), 1, outarray.get(), 1, outarray.get(), 1);
1421 IProductWRTBase_SumFacKernel(base0, base1, dbase2, wsp5, wsp2, wsp0, true,
1422 true, false);
1423 Vmath::Vadd(m_ncoeffs, wsp2.get(), 1, outarray.get(), 1, outarray.get(), 1);
1424}
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.hpp:439

References ASSERTL1, Nektar::LocalRegions::Expansion::ComputeLaplacianMetric(), Nektar::LocalRegions::eMetricLaplacian00, Nektar::LocalRegions::eMetricLaplacian01, Nektar::LocalRegions::eMetricLaplacian02, Nektar::LocalRegions::eMetricLaplacian11, Nektar::LocalRegions::eMetricLaplacian12, Nektar::LocalRegions::eMetricLaplacian22, Nektar::StdRegions::StdExpansion3D::IProductWRTBase_SumFacKernel(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_metrics, Nektar::StdRegions::StdExpansion::m_ncoeffs, Vmath::Vadd(), Vmath::Vvtvp(), and Vmath::Vvtvvtp().

◆ v_MassLevelCurvatureMatrixOp()

void Nektar::LocalRegions::HexExp::v_MassLevelCurvatureMatrixOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdRegions::StdMatrixKey mkey 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1186 of file HexExp.cpp.

1189{
1190 StdExpansion::MassLevelCurvatureMatrixOp_MatFree(inarray, outarray, mkey);
1191}

◆ v_MassMatrixOp()

void Nektar::LocalRegions::HexExp::v_MassMatrixOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdRegions::StdMatrixKey mkey 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdHexExp.

Definition at line 1149 of file HexExp.cpp.

1152{
1153 StdExpansion::MassMatrixOp_MatFree(inarray, outarray, mkey);
1154}

◆ v_NormalTraceDerivFactors()

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

: 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 1468 of file HexExp.cpp.

1472{
1473 int nquad0 = GetNumPoints(0);
1474 int nquad1 = GetNumPoints(1);
1475 int nquad2 = GetNumPoints(2);
1476
1477 const Array<TwoD, const NekDouble> &df =
1478 m_metricinfo->GetDerivFactors(GetPointsKeys());
1479
1480 if (d0factors.size() != 6)
1481 {
1482 d0factors = Array<OneD, Array<OneD, NekDouble>>(6);
1483 d1factors = Array<OneD, Array<OneD, NekDouble>>(6);
1484 d2factors = Array<OneD, Array<OneD, NekDouble>>(6);
1485 }
1486
1487 if (d0factors[0].size() != nquad0 * nquad1)
1488 {
1489 d0factors[0] = Array<OneD, NekDouble>(nquad0 * nquad1);
1490 d0factors[5] = Array<OneD, NekDouble>(nquad0 * nquad1);
1491 d1factors[0] = Array<OneD, NekDouble>(nquad0 * nquad1);
1492 d1factors[5] = Array<OneD, NekDouble>(nquad0 * nquad1);
1493 d2factors[0] = Array<OneD, NekDouble>(nquad0 * nquad1);
1494 d2factors[5] = Array<OneD, NekDouble>(nquad0 * nquad1);
1495 }
1496
1497 if (d0factors[1].size() != nquad0 * nquad2)
1498 {
1499 d0factors[1] = Array<OneD, NekDouble>(nquad0 * nquad2);
1500 d0factors[3] = Array<OneD, NekDouble>(nquad0 * nquad2);
1501 d1factors[1] = Array<OneD, NekDouble>(nquad0 * nquad2);
1502 d1factors[3] = Array<OneD, NekDouble>(nquad0 * nquad2);
1503 d2factors[1] = Array<OneD, NekDouble>(nquad0 * nquad2);
1504 d2factors[3] = Array<OneD, NekDouble>(nquad0 * nquad2);
1505 }
1506
1507 if (d0factors[2].size() != nquad1 * nquad2)
1508 {
1509 d0factors[2] = Array<OneD, NekDouble>(nquad1 * nquad2);
1510 d0factors[4] = Array<OneD, NekDouble>(nquad1 * nquad2);
1511 d1factors[2] = Array<OneD, NekDouble>(nquad1 * nquad2);
1512 d1factors[4] = Array<OneD, NekDouble>(nquad1 * nquad2);
1513 d2factors[2] = Array<OneD, NekDouble>(nquad1 * nquad2);
1514 d2factors[4] = Array<OneD, NekDouble>(nquad1 * nquad2);
1515 }
1516
1517 // Outwards normals
1518 const Array<OneD, const Array<OneD, NekDouble>> &normal_0 =
1519 GetTraceNormal(0);
1520 const Array<OneD, const Array<OneD, NekDouble>> &normal_1 =
1521 GetTraceNormal(1);
1522 const Array<OneD, const Array<OneD, NekDouble>> &normal_2 =
1523 GetTraceNormal(2);
1524 const Array<OneD, const Array<OneD, NekDouble>> &normal_3 =
1525 GetTraceNormal(3);
1526 const Array<OneD, const Array<OneD, NekDouble>> &normal_4 =
1527 GetTraceNormal(4);
1528 const Array<OneD, const Array<OneD, NekDouble>> &normal_5 =
1529 GetTraceNormal(5);
1530
1531 int ncoords = normal_0.size();
1532
1533 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1534 {
1535 // faces 0 and 5
1536 for (int i = 0; i < nquad0 * nquad1; ++i)
1537 {
1538 d0factors[0][i] = df[0][i] * normal_0[0][i];
1539 d1factors[0][i] = df[1][i] * normal_0[0][i];
1540 d2factors[0][i] = df[2][i] * normal_0[0][i];
1541
1542 d0factors[5][i] =
1543 df[0][nquad0 * nquad1 * (nquad2 - 1) + i] * normal_5[0][i];
1544 d1factors[5][i] =
1545 df[1][nquad0 * nquad1 * (nquad2 - 1) + i] * normal_5[0][i];
1546 d2factors[5][i] =
1547 df[2][nquad0 * nquad1 * (nquad2 - 1) + i] * normal_5[0][i];
1548 }
1549
1550 for (int n = 1; n < ncoords; ++n)
1551 {
1552 for (int i = 0; i < nquad0 * nquad1; ++i)
1553 {
1554 d0factors[0][i] += df[3 * n][i] * normal_0[n][i];
1555 d1factors[0][i] += df[3 * n + 1][i] * normal_0[n][i];
1556 d2factors[0][i] += df[3 * n + 2][i] * normal_0[n][i];
1557
1558 d0factors[5][i] +=
1559 df[3 * n][nquad0 * nquad1 * (nquad2 - 1) + i] *
1560 normal_5[n][i];
1561 d1factors[5][i] +=
1562 df[3 * n + 1][nquad0 * nquad1 * (nquad2 - 1) + i] *
1563 normal_5[n][i];
1564 d2factors[5][i] +=
1565 df[3 * n + 2][nquad0 * nquad1 * (nquad2 - 1) + i] *
1566 normal_5[n][i];
1567 }
1568 }
1569
1570 // faces 1 and 3
1571 for (int j = 0; j < nquad2; ++j)
1572 {
1573 for (int i = 0; i < nquad0; ++i)
1574 {
1575 d0factors[1][j * nquad0 + i] = df[0][j * nquad0 * nquad1 + i] *
1576 normal_1[0][j * nquad0 + i];
1577 d1factors[1][j * nquad0 + i] = df[1][j * nquad0 * nquad1 + i] *
1578 normal_1[0][j * nquad0 + i];
1579 d2factors[1][j * nquad0 + i] = df[2][j * nquad0 * nquad1 + i] *
1580 normal_1[0][j * nquad0 + i];
1581
1582 d0factors[3][j * nquad0 + i] =
1583 df[0][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1584 normal_3[0][j * nquad0 + i];
1585 d1factors[3][j * nquad0 + i] =
1586 df[1][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1587 normal_3[0][j * nquad0 + i];
1588 d2factors[3][j * nquad0 + i] =
1589 df[2][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1590 normal_3[0][j * nquad0 + i];
1591 }
1592 }
1593
1594 for (int n = 1; n < ncoords; ++n)
1595 {
1596 for (int j = 0; j < nquad2; ++j)
1597 {
1598 for (int i = 0; i < nquad0; ++i)
1599 {
1600 d0factors[1][j * nquad0 + i] +=
1601 df[3 * n][j * nquad0 * nquad1 + i] *
1602 normal_1[0][j * nquad0 + i];
1603 d1factors[1][j * nquad0 + i] +=
1604 df[3 * n + 1][j * nquad0 * nquad1 + i] *
1605 normal_1[0][j * nquad0 + i];
1606 d2factors[1][j * nquad0 + i] +=
1607 df[3 * n + 2][j * nquad0 * nquad1 + i] *
1608 normal_1[0][j * nquad0 + i];
1609
1610 d0factors[3][j * nquad0 + i] +=
1611 df[3 * n][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1612 normal_3[0][j * nquad0 + i];
1613 d1factors[3][j * nquad0 + i] +=
1614 df[3 * n + 1][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1615 normal_3[0][j * nquad0 + i];
1616 d2factors[3][j * nquad0 + i] +=
1617 df[3 * n + 2][(j + 1) * nquad0 * nquad1 - nquad0 + i] *
1618 normal_3[0][j * nquad0 + i];
1619 }
1620 }
1621 }
1622
1623 // faces 2 and 4
1624 for (int j = 0; j < nquad2; ++j)
1625 {
1626 for (int i = 0; i < nquad1; ++i)
1627 {
1628 d0factors[2][j * nquad1 + i] =
1629 df[0][j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1630 normal_2[0][j * nquad1 + i];
1631 d1factors[2][j * nquad1 + i] =
1632 df[1][j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1633 normal_2[0][j * nquad1 + i];
1634 d2factors[2][j * nquad1 + i] =
1635 df[2][j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1636 normal_2[0][j * nquad1 + i];
1637
1638 d0factors[4][j * nquad1 + i] =
1639 df[0][j * nquad0 * nquad1 + i * nquad0] *
1640 normal_4[0][j * nquad1 + i];
1641 d1factors[4][j * nquad1 + i] =
1642 df[1][j * nquad0 * nquad1 + i * nquad0] *
1643 normal_4[0][j * nquad1 + i];
1644 d2factors[4][j * nquad1 + i] =
1645 df[2][j * nquad0 * nquad1 + i * nquad0] *
1646 normal_4[0][j * nquad1 + i];
1647 }
1648 }
1649
1650 for (int n = 1; n < ncoords; ++n)
1651 {
1652 for (int j = 0; j < nquad2; ++j)
1653 {
1654 for (int i = 0; i < nquad1; ++i)
1655 {
1656 d0factors[2][j * nquad1 + i] +=
1657 df[3 * n][j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1658 normal_2[n][j * nquad1 + i];
1659 d1factors[2][j * nquad1 + i] +=
1660 df[3 * n + 1]
1661 [j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1662 normal_2[n][j * nquad1 + i];
1663 d2factors[2][j * nquad1 + i] +=
1664 df[3 * n + 2]
1665 [j * nquad0 * nquad1 + (i + 1) * nquad0 - 1] *
1666 normal_2[n][j * nquad1 + i];
1667
1668 d0factors[4][j * nquad1 + i] +=
1669 df[3 * n][i * nquad0 + j * nquad0 * nquad1] *
1670 normal_4[n][j * nquad1 + i];
1671 d1factors[4][j * nquad1 + i] +=
1672 df[3 * n + 1][i * nquad0 + j * nquad0 * nquad1] *
1673 normal_4[n][j * nquad1 + i];
1674 d2factors[4][j * nquad1 + i] +=
1675 df[3 * n + 2][i * nquad0 + j * nquad0 * nquad1] *
1676 normal_4[n][j * nquad1 + i];
1677 }
1678 }
1679 }
1680 }
1681 else
1682 {
1683 // Faces 0 and 5
1684 for (int i = 0; i < nquad0 * nquad1; ++i)
1685 {
1686 d0factors[0][i] = df[0][0] * normal_0[0][i];
1687 d0factors[5][i] = df[0][0] * normal_5[0][i];
1688
1689 d1factors[0][i] = df[1][0] * normal_0[0][i];
1690 d1factors[5][i] = df[1][0] * normal_5[0][i];
1691
1692 d2factors[0][i] = df[2][0] * normal_0[0][i];
1693 d2factors[5][i] = df[2][0] * normal_5[0][i];
1694 }
1695
1696 for (int n = 1; n < ncoords; ++n)
1697 {
1698 for (int i = 0; i < nquad0 * nquad1; ++i)
1699 {
1700 d0factors[0][i] += df[3 * n][0] * normal_0[n][i];
1701 d0factors[5][i] += df[3 * n][0] * normal_5[n][i];
1702
1703 d1factors[0][i] += df[3 * n + 1][0] * normal_0[n][i];
1704 d1factors[5][i] += df[3 * n + 1][0] * normal_5[n][i];
1705
1706 d2factors[0][i] += df[3 * n + 2][0] * normal_0[n][i];
1707 d2factors[5][i] += df[3 * n + 2][0] * normal_5[n][i];
1708 }
1709 }
1710
1711 // faces 1 and 3
1712 for (int i = 0; i < nquad0 * nquad2; ++i)
1713 {
1714 d0factors[1][i] = df[0][0] * normal_1[0][i];
1715 d0factors[3][i] = df[0][0] * normal_3[0][i];
1716
1717 d1factors[1][i] = df[1][0] * normal_1[0][i];
1718 d1factors[3][i] = df[1][0] * normal_3[0][i];
1719
1720 d2factors[1][i] = df[2][0] * normal_1[0][i];
1721 d2factors[3][i] = df[2][0] * normal_3[0][i];
1722 }
1723
1724 for (int n = 1; n < ncoords; ++n)
1725 {
1726 for (int i = 0; i < nquad0 * nquad2; ++i)
1727 {
1728 d0factors[1][i] += df[3 * n][0] * normal_1[n][i];
1729 d0factors[3][i] += df[3 * n][0] * normal_3[n][i];
1730
1731 d1factors[1][i] += df[3 * n + 1][0] * normal_1[n][i];
1732 d1factors[3][i] += df[3 * n + 1][0] * normal_3[n][i];
1733
1734 d2factors[1][i] += df[3 * n + 2][0] * normal_1[n][i];
1735 d2factors[3][i] += df[3 * n + 2][0] * normal_3[n][i];
1736 }
1737 }
1738
1739 // faces 2 and 4
1740 for (int i = 0; i < nquad1 * nquad2; ++i)
1741 {
1742 d0factors[2][i] = df[0][0] * normal_2[0][i];
1743 d0factors[4][i] = df[0][0] * normal_4[0][i];
1744
1745 d1factors[2][i] = df[1][0] * normal_2[0][i];
1746 d1factors[4][i] = df[1][0] * normal_4[0][i];
1747
1748 d2factors[2][i] = df[2][0] * normal_2[0][i];
1749 d2factors[4][i] = df[2][0] * normal_4[0][i];
1750 }
1751
1752 for (int n = 1; n < ncoords; ++n)
1753 {
1754 for (int i = 0; i < nquad1 * nquad2; ++i)
1755 {
1756 d0factors[2][i] += df[3 * n][0] * normal_2[n][i];
1757 d0factors[4][i] += df[3 * n][0] * normal_4[n][i];
1758
1759 d1factors[2][i] += df[3 * n + 1][0] * normal_2[n][i];
1760 d1factors[4][i] += df[3 * n + 1][0] * normal_4[n][i];
1761
1762 d2factors[2][i] += df[3 * n + 2][0] * normal_2[n][i];
1763 d2factors[4][i] += df[3 * n + 2][0] * normal_4[n][i];
1764 }
1765 }
1766 }
1767}
const NormalVector & GetTraceNormal(const int id)
Definition: Expansion.cpp:251

References Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetNumPoints(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::LocalRegions::Expansion::GetTraceNormal(), and Nektar::LocalRegions::Expansion::m_metricinfo.

◆ v_PhysDeriv() [1/2]

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

Calculate the derivative of the physical points.

For Hexahedral region can use the Tensor_Deriv function defined under StdExpansion.

Parameters
inarrayInput array
out_d0Derivative of inarray in first direction.
out_d1Derivative of inarray in second direction.
out_d2Derivative of inarray in third direction.

Reimplemented from Nektar::StdRegions::StdHexExp.

Definition at line 142 of file HexExp.cpp.

146{
147 int nquad0 = m_base[0]->GetNumPoints();
148 int nquad1 = m_base[1]->GetNumPoints();
149 int nquad2 = m_base[2]->GetNumPoints();
150 int ntot = nquad0 * nquad1 * nquad2;
151
152 Array<TwoD, const NekDouble> df =
153 m_metricinfo->GetDerivFactors(GetPointsKeys());
154 Array<OneD, NekDouble> Diff0 = Array<OneD, NekDouble>(ntot);
155 Array<OneD, NekDouble> Diff1 = Array<OneD, NekDouble>(ntot);
156 Array<OneD, NekDouble> Diff2 = Array<OneD, NekDouble>(ntot);
157
158 StdHexExp::v_PhysDeriv(inarray, Diff0, Diff1, Diff2);
159
160 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
161 {
162 if (out_d0.size())
163 {
164 Vmath::Vmul(ntot, &df[0][0], 1, &Diff0[0], 1, &out_d0[0], 1);
165 Vmath::Vvtvp(ntot, &df[1][0], 1, &Diff1[0], 1, &out_d0[0], 1,
166 &out_d0[0], 1);
167 Vmath::Vvtvp(ntot, &df[2][0], 1, &Diff2[0], 1, &out_d0[0], 1,
168 &out_d0[0], 1);
169 }
170
171 if (out_d1.size())
172 {
173 Vmath::Vmul(ntot, &df[3][0], 1, &Diff0[0], 1, &out_d1[0], 1);
174 Vmath::Vvtvp(ntot, &df[4][0], 1, &Diff1[0], 1, &out_d1[0], 1,
175 &out_d1[0], 1);
176 Vmath::Vvtvp(ntot, &df[5][0], 1, &Diff2[0], 1, &out_d1[0], 1,
177 &out_d1[0], 1);
178 }
179
180 if (out_d2.size())
181 {
182 Vmath::Vmul(ntot, &df[6][0], 1, &Diff0[0], 1, &out_d2[0], 1);
183 Vmath::Vvtvp(ntot, &df[7][0], 1, &Diff1[0], 1, &out_d2[0], 1,
184 &out_d2[0], 1);
185 Vmath::Vvtvp(ntot, &df[8][0], 1, &Diff2[0], 1, &out_d2[0], 1,
186 &out_d2[0], 1);
187 }
188 }
189 else // regular geometry
190 {
191 if (out_d0.size())
192 {
193 Vmath::Smul(ntot, df[0][0], &Diff0[0], 1, &out_d0[0], 1);
194 Blas::Daxpy(ntot, df[1][0], &Diff1[0], 1, &out_d0[0], 1);
195 Blas::Daxpy(ntot, df[2][0], &Diff2[0], 1, &out_d0[0], 1);
196 }
197
198 if (out_d1.size())
199 {
200 Vmath::Smul(ntot, df[3][0], &Diff0[0], 1, &out_d1[0], 1);
201 Blas::Daxpy(ntot, df[4][0], &Diff1[0], 1, &out_d1[0], 1);
202 Blas::Daxpy(ntot, df[5][0], &Diff2[0], 1, &out_d1[0], 1);
203 }
204
205 if (out_d2.size())
206 {
207 Vmath::Smul(ntot, df[6][0], &Diff0[0], 1, &out_d2[0], 1);
208 Blas::Daxpy(ntot, df[7][0], &Diff1[0], 1, &out_d2[0], 1);
209 Blas::Daxpy(ntot, df[8][0], &Diff2[0], 1, &out_d2[0], 1);
210 }
211 }
212}
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:135

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_PhysDeriv() [2/2]

void Nektar::LocalRegions::HexExp::v_PhysDeriv ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Calculate the derivative of the physical points in a single direction.

Parameters
dirDirection in which to compute derivative. Valid values are 0, 1, 2.
inarrayInput array.
outarrayOutput array.

Reimplemented from Nektar::StdRegions::StdHexExp.

Definition at line 223 of file HexExp.cpp.

226{
227 switch (dir)
228 {
229 case 0:
230 {
231 PhysDeriv(inarray, outarray, NullNekDouble1DArray,
233 }
234 break;
235 case 1:
236 {
237 PhysDeriv(inarray, NullNekDouble1DArray, outarray,
239 }
240 break;
241 case 2:
242 {
244 outarray);
245 }
246 break;
247 default:
248 {
249 ASSERTL1(false, "input dir is out of range");
250 }
251 break;
252 }
253}
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)
Definition: StdExpansion.h:849
static Array< OneD, NekDouble > NullNekDouble1DArray

References ASSERTL1, Nektar::NullNekDouble1DArray, and Nektar::StdRegions::StdExpansion::PhysDeriv().

◆ v_PhysDirectionalDeriv()

void Nektar::LocalRegions::HexExp::v_PhysDirectionalDeriv ( const Array< OneD, const NekDouble > &  inarray,
const Array< OneD, const NekDouble > &  direction,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Physical derivative along a direction vector.

See also
StdRegions::StdExpansion::PhysDirectionalDeriv

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 255 of file HexExp.cpp.

259{
260
261 int shapedim = 3;
262 int nquad0 = m_base[0]->GetNumPoints();
263 int nquad1 = m_base[1]->GetNumPoints();
264 int nquad2 = m_base[2]->GetNumPoints();
265 int ntot = nquad0 * nquad1 * nquad2;
266
267 Array<TwoD, const NekDouble> df =
268 m_metricinfo->GetDerivFactors(GetPointsKeys());
269 Array<OneD, NekDouble> Diff0 = Array<OneD, NekDouble>(ntot);
270 Array<OneD, NekDouble> Diff1 = Array<OneD, NekDouble>(ntot);
271 Array<OneD, NekDouble> Diff2 = Array<OneD, NekDouble>(ntot);
272
273 StdHexExp::v_PhysDeriv(inarray, Diff0, Diff1, Diff2);
274
275 Array<OneD, Array<OneD, NekDouble>> dfdir(shapedim);
276 Expansion::ComputeGmatcdotMF(df, direction, dfdir);
277
278 Vmath::Vmul(ntot, &dfdir[0][0], 1, &Diff0[0], 1, &outarray[0], 1);
279 Vmath::Vvtvp(ntot, &dfdir[1][0], 1, &Diff1[0], 1, &outarray[0], 1,
280 &outarray[0], 1);
281 Vmath::Vvtvp(ntot, &dfdir[2][0], 1, &Diff2[0], 1, &outarray[0], 1,
282 &outarray[0], 1);
283}

References Nektar::LocalRegions::Expansion::ComputeGmatcdotMF(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_metricinfo, Vmath::Vmul(), and Vmath::Vvtvp().

◆ v_PhysEvaluate() [1/2]

NekDouble Nektar::LocalRegions::HexExp::v_PhysEvaluate ( const Array< OneD, const NekDouble > &  coords,
const Array< OneD, const NekDouble > &  physvals 
)
overrideprotectedvirtual

This function evaluates the expansion at a single (arbitrary) point of the domain.

Based on the value of the expansion at the quadrature points, this function calculates the value of the expansion at an arbitrary single points (with coordinates \( \mathbf{x_c}\) given by the pointer coords). This operation, equivalent to

\[ u(\mathbf{x_c}) = \sum_p \phi_p(\mathbf{x_c}) \hat{u}_p \]

is evaluated using Lagrangian interpolants through the quadrature points:

\[ u(\mathbf{x_c}) = \sum_p h_p(\mathbf{x_c}) u_p\]

This function requires that the physical value array \(\mathbf{u}\) (implemented as the attribute #phys) is set.

Parameters
coordsthe coordinates of the single point
Returns
returns the value of the expansion at the single point

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 581 of file HexExp.cpp.

583{
584 Array<OneD, NekDouble> Lcoord = Array<OneD, NekDouble>(3);
585
586 ASSERTL0(m_geom, "m_geom not defined");
587 m_geom->GetLocCoords(coord, Lcoord);
588 return StdExpansion3D::v_PhysEvaluate(Lcoord, physvals);
589}

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

◆ v_PhysEvaluate() [2/2]

NekDouble Nektar::LocalRegions::HexExp::v_PhysEvaluate ( const Array< OneD, NekDouble > &  coord,
const Array< OneD, const NekDouble > &  inarray,
std::array< NekDouble, 3 > &  firstOrderDerivs 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdHexExp.

Definition at line 591 of file HexExp.cpp.

594{
595 Array<OneD, NekDouble> Lcoord(3);
596 ASSERTL0(m_geom, "m_geom not defined");
597 m_geom->GetLocCoords(coord, Lcoord);
598 return StdHexExp::v_PhysEvaluate(Lcoord, inarray, firstOrderDerivs);
599}

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

◆ v_ReduceOrderCoeffs()

void Nektar::LocalRegions::HexExp::v_ReduceOrderCoeffs ( int  numMin,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

This function is used to compute exactly the advective numerical flux on the interface of two elements with different expansions, hence an appropriate number of Gauss points has to be used. The number of Gauss points has to be equal to the number used by the highest polynomial degree of the two adjacent elements

Parameters
numMinIs the reduced polynomial order
inarrayInput array of coefficients
dumpVarOutput array of reduced coefficients.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1211 of file HexExp.cpp.

1214{
1215 int n_coeffs = inarray.size();
1216 int nmodes0 = m_base[0]->GetNumModes();
1217 int nmodes1 = m_base[1]->GetNumModes();
1218 int nmodes2 = m_base[2]->GetNumModes();
1219 int numMax = nmodes0;
1220
1221 Array<OneD, NekDouble> coeff(n_coeffs);
1222 Array<OneD, NekDouble> coeff_tmp1(nmodes0 * nmodes1, 0.0);
1223 Array<OneD, NekDouble> coeff_tmp2(n_coeffs, 0.0);
1224 Array<OneD, NekDouble> tmp, tmp2, tmp3, tmp4;
1225
1226 Vmath::Vcopy(n_coeffs, inarray, 1, coeff_tmp2, 1);
1227
1228 const LibUtilities::PointsKey Pkey0(nmodes0,
1230 const LibUtilities::PointsKey Pkey1(nmodes1,
1232 const LibUtilities::PointsKey Pkey2(nmodes2,
1234
1235 LibUtilities::BasisKey b0(m_base[0]->GetBasisType(), nmodes0, Pkey0);
1236 LibUtilities::BasisKey b1(m_base[1]->GetBasisType(), nmodes1, Pkey1);
1237 LibUtilities::BasisKey b2(m_base[2]->GetBasisType(), nmodes2, Pkey2);
1238 LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A, nmodes0, Pkey0);
1239 LibUtilities::BasisKey bortho1(LibUtilities::eOrtho_A, nmodes1, Pkey1);
1240 LibUtilities::BasisKey bortho2(LibUtilities::eOrtho_A, nmodes2, Pkey2);
1241
1242 LibUtilities::InterpCoeff3D(b0, b1, b2, coeff_tmp2, bortho0, bortho1,
1243 bortho2, coeff);
1244
1245 Vmath::Zero(n_coeffs, coeff_tmp2, 1);
1246
1247 int cnt = 0, cnt2 = 0;
1248
1249 for (int u = 0; u < numMin + 1; ++u)
1250 {
1251 for (int i = 0; i < numMin; ++i)
1252 {
1253 Vmath::Vcopy(numMin, tmp = coeff + cnt + cnt2, 1,
1254 tmp2 = coeff_tmp1 + cnt, 1);
1255
1256 cnt = i * numMax;
1257 }
1258
1259 Vmath::Vcopy(nmodes0 * nmodes1, tmp3 = coeff_tmp1, 1,
1260 tmp4 = coeff_tmp2 + cnt2, 1);
1261
1262 cnt2 = u * nmodes0 * nmodes1;
1263 }
1264
1265 LibUtilities::InterpCoeff3D(bortho0, bortho1, bortho2, coeff_tmp2, b0, b1,
1266 b2, outarray);
1267}
void InterpCoeff3D(const BasisKey &fbasis0, const BasisKey &fbasis1, const BasisKey &fbasis2, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, const BasisKey &tbasis2, Array< OneD, NekDouble > &to)
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:42

References Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eOrtho_A, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::LibUtilities::InterpCoeff3D(), Nektar::StdRegions::StdExpansion::m_base, Vmath::Vcopy(), and Vmath::Zero().

◆ v_StdPhysEvaluate()

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

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 573 of file HexExp.cpp.

576{
577 // Evaluate point in local coordinates.
578 return StdExpansion3D::v_PhysEvaluate(Lcoord, physvals);
579}

◆ v_SVVLaplacianFilter()

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

Reimplemented from Nektar::StdRegions::StdHexExp.

Definition at line 1269 of file HexExp.cpp.

1271{
1272 int nq = GetTotPoints();
1273
1274 // Calculate sqrt of the Jacobian
1275 Array<OneD, const NekDouble> jac = m_metricinfo->GetJac(GetPointsKeys());
1276 Array<OneD, NekDouble> sqrt_jac(nq);
1277 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1278 {
1279 Vmath::Vsqrt(nq, jac, 1, sqrt_jac, 1);
1280 }
1281 else
1282 {
1283 Vmath::Fill(nq, sqrt(jac[0]), sqrt_jac, 1);
1284 }
1285
1286 // Multiply array by sqrt(Jac)
1287 Vmath::Vmul(nq, sqrt_jac, 1, array, 1, array, 1);
1288
1289 // Apply std region filter
1290 StdHexExp::v_SVVLaplacianFilter(array, mkey);
1291
1292 // Divide by sqrt(Jac)
1293 Vmath::Vdiv(nq, array, 1, sqrt_jac, 1, array, 1);
1294}
void Vdiv(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x/y.
Definition: Vmath.hpp:126

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

◆ v_WeakDerivMatrixOp()

void Nektar::LocalRegions::HexExp::v_WeakDerivMatrixOp ( const int  i,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdRegions::StdMatrixKey mkey 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdHexExp.

Definition at line 1171 of file HexExp.cpp.

1175{
1176 StdExpansion::WeakDerivMatrixOp_MatFree(i, inarray, outarray, mkey);
1177}

◆ v_WeakDirectionalDerivMatrixOp()

void Nektar::LocalRegions::HexExp::v_WeakDirectionalDerivMatrixOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdRegions::StdMatrixKey mkey 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1179 of file HexExp.cpp.

1182{
1183 StdExpansion::WeakDirectionalDerivMatrixOp_MatFree(inarray, outarray, mkey);
1184}

Member Data Documentation

◆ m_matrixManager

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

Definition at line 244 of file HexExp.h.

Referenced by v_DropLocMatrix(), v_FwdTrans(), and v_GetLocMatrix().

◆ m_staticCondMatrixManager

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

Definition at line 246 of file HexExp.h.

Referenced by v_DropLocStaticCondMatrix(), and v_GetLocStaticCondMatrix().