Nektar++
Loading...
Searching...
No Matches
Public Member Functions | Protected Member Functions | List of all members
Nektar::StdRegions::StdHexExp Class Reference

Class representing a hexehedral element in reference space. More...

#include <StdHexExp.h>

Inheritance diagram for Nektar::StdRegions::StdHexExp:
[legend]

Public Member Functions

 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.
 
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
 
int GetEdgeNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th edge.
 
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.
 
 StdExpansion (const int numcoeffs, const int numbases, const LibUtilities::BasisKey &Ba=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bb=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bc=LibUtilities::NullBasisKey)
 Constructor.
 
 StdExpansion (const StdExpansion &T)
 Copy Constructor.
 
virtual ~StdExpansion ()
 Destructor.
 
int GetNumBases () const
 This function returns the number of 1D bases used in the expansion.
 
const Array< OneD, const LibUtilities::BasisSharedPtr > & GetBase () const
 This function gets the shared point to basis.
 
const LibUtilities::BasisSharedPtrGetBasis (int dir) const
 This function gets the shared point to basis in the dir direction.
 
int GetNcoeffs (void) const
 This function returns the total number of coefficients used in the expansion.
 
int GetTotPoints () const
 This function returns the total number of quadrature points used in the element.
 
LibUtilities::BasisType GetBasisType (const int dir) const
 This function returns the type of basis used in the dir direction.
 
int GetBasisNumModes (const int dir) const
 This function returns the number of expansion modes in the dir direction.
 
int EvalBasisNumModesMax (void) const
 This function returns the maximum number of expansion modes over all local directions.
 
LibUtilities::PointsType GetPointsType (const int dir) const
 This function returns the type of quadrature points used in the dir direction.
 
int GetNumPoints (const int dir) const
 This function returns the number of quadrature points in the dir direction.
 
const Array< OneD, const NekDouble > & GetPoints (const int dir) const
 This function returns a pointer to the array containing the quadrature points in dir direction.
 
int GetNverts () const
 This function returns the number of vertices of the expansion domain.
 
int GetTraceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th trace.
 
int GetTraceIntNcoeffs (const int i) const
 
int GetTraceNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th trace.
 
const LibUtilities::BasisKey GetTraceBasisKey (const int i, int k=-1, bool UseGLL=false) const
 This function returns the basis key belonging to the i-th trace.
 
LibUtilities::PointsKey GetTracePointsKey (const int i, int k=-1) const
 This function returns the basis key belonging to the i-th trace.
 
int NumBndryCoeffs (void) const
 
int NumDGBndryCoeffs (void) const
 
const LibUtilities::PointsKey GetNodalPointsKey () const
 This function returns the type of expansion Nodal point type if defined.
 
int GetNtraces () const
 Returns the number of trace elements connected to this element.
 
LibUtilities::ShapeType DetShapeType () const
 This function returns the shape of the expansion domain.
 
std::shared_ptr< StdExpansionGetStdExp () const
 
std::shared_ptr< StdExpansionGetLinStdExp (void) const
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion () const
 
bool IsNodalNonTensorialExp ()
 
void NodalToModal (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space.
 
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space.
 
void FwdTransBndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble Integral (const Array< OneD, const NekDouble > &inarray)
 This function integrates the specified function over the domain.
 
void FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 This function fills the array outarray with the mode-th mode of the expansion.
 
void IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 this function calculates the inner product of a given function f with the different modes of the expansion
 
void IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check)
 
void IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
int GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id.
 
void SetElmtId (const int id)
 Set the element id of this expansion when used in a list by returning value of m_elmt_id.
 
void GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2=NullNekDouble1DArray, Array< OneD, NekDouble > &coords_3=NullNekDouble1DArray)
 this function returns the physical coordinates of the quadrature points of the expansion
 
void GetCoord (const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
 given the coordinates of a point of the element in the local collapsed coordinate system, this function calculates the physical coordinates of the point
 
DNekMatSharedPtr GetStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr GetStdStaticCondMatrix (const StdMatrixKey &mkey)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
void DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
int CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
NekDouble StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
int GetCoordim ()
 
void GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
void GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
int GetVertexMap (const int localVertexId, bool useCoeffPacking=false)
 
void GetTraceToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
 
void GetTraceCoeffMap (const unsigned int traceid, Array< OneD, unsigned int > &maparray)
 
void GetElmtTraceToTraceMap (const unsigned int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards, int P=-1, int Q=-1)
 
void GetTraceInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, const Orientation traceOrient=eForwards)
 
void GetTraceNumModes (const int tid, int &numModes0, int &numModes1, const Orientation traceOrient=eDir1FwdDir1_Dir2FwdDir2)
 
void MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
DNekMatSharedPtr CreateGeneralMatrix (const StdMatrixKey &mkey)
 this function generates the mass matrix \(\mathbf{M}[i][j] = \int \phi_i(\mathbf{x}) \phi_j(\mathbf{x}) d\mathbf{x}\)
 
void GeneralMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
void ExponentialFilter (Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
 
void LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
DNekMatSharedPtr GenMatrix (const StdMatrixKey &mkey)
 
void PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void PhysDeriv_s (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds)
 
void PhysDeriv_n (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn)
 
void PhysDirectionalDeriv (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &outarray)
 
void StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
NekDouble PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain.
 
NekDouble PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
 This function evaluates the first derivative of the expansion at a single (arbitrary) point of the domain.
 
NekDouble PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs, std::array< NekDouble, 6 > &secondOrderDerivs)
 
NekDouble PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain.
 
NekDouble PhysEvaluateBasis (const Array< OneD, const NekDouble > &coords, int mode)
 This function evaluates the basis function mode mode at a point coords of the domain.
 
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta.
 
void LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi)
 Convert local collapsed coordinates eta into local cartesian coordinate xi.
 
void PhysInterp (std::shared_ptr< StdExpansion > fromExp, const Array< OneD, const NekDouble > &fromData, Array< OneD, NekDouble > &toData)
 interpolate from one set of quadrature points available from FromExp to the set of quadrature points in the current expansion. If the points are the same this routine will just copy the data
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual void v_DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
NekDouble Linf (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( L_\infty\) error \( |\epsilon|_\infty = \max |u - u_{exact}|\) where \( u_{exact}\) is given by the array sol.
 
NekDouble L2 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( L_2\) error, \( | \epsilon |_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 dx \right]^{1/2} d\xi_1 \) where \( u_{exact}\) is given by the array sol.
 
NekDouble H1 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete \( H^1\) error, \( | \epsilon |^1_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 + \nabla(u - u_{exact})\cdot\nabla(u - u_{exact})\cdot dx \right]^{1/2} d\xi_1 \) where \( u_{exact}\) is given by the array sol.
 
const LibUtilities::PointsKeyVector GetPointsKeys () const
 
DNekMatSharedPtr BuildInverseTransformationMatrix (const DNekScalMatSharedPtr &m_transformationmatrix)
 
void PhysInterpToSimplexEquiSpaced (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int npset=-1)
 This function performs an interpolation from the physical space points provided at input into an array of equispaced points which are not the collapsed coordinate. So for a tetrahedron you will only get a tetrahedral number of values.
 
void GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true)
 This function provides the connectivity of local simplices (triangles or tets) to connect the equispaced data points provided by PhysInterpToSimplexEquiSpaced.
 
void EquiSpacedToCoeffs (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs a projection/interpolation from the equispaced points sometimes used in post-processing onto the coefficient space.
 
template<class T >
std::shared_ptr< T > as ()
 
void IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
void GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 

Protected Member Functions

void 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.
 
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_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_PhysEvalFirstDeriv (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, bool useGLL=false) 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.
 
NekDouble v_PhysEvaluateInterp (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals) override
 
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.
 
NekDouble BaryTensorDeriv (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
 
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
 
void v_PhysInterp (std::shared_ptr< StdExpansion > fromExp, const Array< OneD, const NekDouble > &fromData, Array< OneD, NekDouble > &toData) 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.
 
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 NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals, NekDouble &deriv, NekDouble &deriv2)
 This function performs the barycentric interpolation of the polynomial stored in coord at a point physvals using barycentric interpolation weights in direction.
 
template<int DIR>
NekDouble BaryEvaluateBasis (const NekDouble &coord, const int &mode)
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals)
 Helper function to pass an unused value by reference into BaryEvaluate.
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals, NekDouble &deriv)
 

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
 

Detailed Description

Class representing a hexehedral element in reference space.

Definition at line 43 of file StdHexExp.h.

Constructor & Destructor Documentation

◆ StdHexExp() [1/2]

Nektar::StdRegions::StdHexExp::StdHexExp ( const LibUtilities::BasisKey Ba,
const LibUtilities::BasisKey Bb,
const LibUtilities::BasisKey Bc 
)

Definition at line 46 of file StdHexExp.cpp.

49 : StdExpansion(Ba.GetNumModes() * Bb.GetNumModes() * Bc.GetNumModes(), 3,
50 Ba, Bb, Bc),
51 StdExpansion3D(Ba.GetNumModes() * Bb.GetNumModes() * Bc.GetNumModes(), Ba,
52 Bb, Bc)
53{
54}
StdExpansion()
Default Constructor.

◆ StdHexExp() [2/2]

Nektar::StdRegions::StdHexExp::StdHexExp ( const StdHexExp T)
default

◆ ~StdHexExp()

Nektar::StdRegions::StdHexExp::~StdHexExp ( )
overridedefault

Member Function Documentation

◆ v_BwdTrans()

void Nektar::StdRegions::StdHexExp::v_BwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Backward transformation is three dimensional tensorial expansion \( u (\xi_{1i}, \xi_{2j}, \xi_{3k}) = \sum_{p=0}^{Q_x} \psi_p^a (\xi_{1i}) \lbrace { \sum_{q=0}^{Q_y} \psi_{q}^a (\xi_{2j}) \lbrace { \sum_{r=0}^{Q_z} \hat u_{pqr} \psi_{r}^a (\xi_{3k}) \rbrace} \rbrace}. \) And sumfactorizing step of the form is as:\ \( f_{r} (\xi_{3k}) = \sum_{r=0}^{Q_z} \hat u_{pqr} \psi_{r}^a (\xi_{3k}),\\ g_{p} (\xi_{2j}, \xi_{3k}) = \sum_{r=0}^{Q_y} \psi_{p}^a (\xi_{2j}) f_{r} (\xi_{3k}),\\ u(\xi_{1i}, \xi_{2j}, \xi_{3k}) = \sum_{p=0}^{Q_x} \psi_{p}^a (\xi_{1i}) g_{p} (\xi_{2j}, \xi_{3k}). \)

Parameters
inarray?
outarray?

Implements Nektar::StdRegions::StdExpansion.

Definition at line 147 of file StdHexExp.cpp.

149{
152 "Basis[1] is not a general tensor type");
153
156 "Basis[2] is not a general tensor type");
157
158 if (m_base[0]->Collocation() && m_base[1]->Collocation() &&
159 m_base[2]->Collocation())
160 {
162 m_base[2]->GetNumPoints(),
163 inarray, 1, outarray, 1);
164 }
165 else
166 {
167 StdHexExp::BwdTrans_SumFac(inarray, outarray);
168 }
169}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
void BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Array< OneD, LibUtilities::BasisSharedPtr > m_base
@ eModified_B
Principle Modified Functions .
Definition BasisType.h:49
@ eModified_C
Principle Modified Functions .
Definition BasisType.h:50
@ eOrtho_C
Principle Orthogonal Functions .
Definition BasisType.h:46
@ eOrtho_B
Principle Orthogonal Functions .
Definition BasisType.h:44
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition Vmath.hpp:825

References ASSERTL1, Nektar::StdRegions::StdExpansion::BwdTrans_SumFac(), Nektar::LibUtilities::eModified_B, Nektar::LibUtilities::eModified_C, Nektar::LibUtilities::eOrtho_B, Nektar::LibUtilities::eOrtho_C, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::GetNumPoints(), Nektar::StdRegions::StdExpansion::m_base, and Vmath::Vcopy().

◆ v_BwdTrans_SumFac()

void Nektar::StdRegions::StdHexExp::v_BwdTrans_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 174 of file StdHexExp.cpp.

176{
177 Array<OneD, NekDouble> wsp(
178 m_base[0]->GetNumPoints() * m_base[2]->GetNumModes() *
179 (m_base[1]->GetNumModes() + m_base[1]->GetNumPoints())); // FIX THIS
180
181 BwdTrans_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetBdata(),
182 m_base[2]->GetBdata(), inarray, outarray, wsp, true,
183 true, true);
184}
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)

References Nektar::StdRegions::StdExpansion3D::BwdTrans_SumFacKernel(), Nektar::StdRegions::StdExpansion::GetNumPoints(), and Nektar::StdRegions::StdExpansion::m_base.

◆ v_BwdTrans_SumFacKernel()

void Nektar::StdRegions::StdHexExp::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 
)
overrideprotectedvirtual
Parameters
base0x-dirn basis matrix
base1y-dirn basis matrix
base2z-dirn basis matrix
inarrayInput vector of modes.
outarrayOutput vector of physical space data.
wspWorkspace of size Q_x*P_z*(P_y+Q_y)
doCheckCollDir0Check for collocation of basis.
doCheckCollDir1Check for collocation of basis.
doCheckCollDir2Check for collocation of basis.
Todo:
Account for some directions being collocated. See StdQuadExp as an example.

Implements Nektar::StdRegions::StdExpansion3D.

Definition at line 199 of file StdHexExp.cpp.

206{
207 int nquad0 = m_base[0]->GetNumPoints();
208 int nquad1 = m_base[1]->GetNumPoints();
209 int nquad2 = m_base[2]->GetNumPoints();
210 int nmodes0 = m_base[0]->GetNumModes();
211 int nmodes1 = m_base[1]->GetNumModes();
212 int nmodes2 = m_base[2]->GetNumModes();
213
214 // Check if using collocation, if requested.
215 bool colldir0 = doCheckCollDir0 ? (m_base[0]->Collocation()) : false;
216 bool colldir1 = doCheckCollDir1 ? (m_base[1]->Collocation()) : false;
217 bool colldir2 = doCheckCollDir2 ? (m_base[2]->Collocation()) : false;
218
219 // If collocation in all directions, Physical values at quadrature
220 // points is just a copy of the modes.
221 if (colldir0 && colldir1 && colldir2)
222 {
223 Vmath::Vcopy(m_ncoeffs, inarray.data(), 1, outarray.data(), 1);
224 }
225 else
226 {
227 // Check sufficiently large workspace.
228 ASSERTL1(wsp.size() >= nquad0 * nmodes2 * (nmodes1 + nquad1),
229 "Workspace size is not sufficient");
230
231 // Assign second half of workspace for 2nd DGEMM operation.
232 Array<OneD, NekDouble> wsp2 = wsp + nquad0 * nmodes1 * nmodes2;
233
234 // BwdTrans in each direction using DGEMM
235 Blas::Dgemm('T', 'T', nmodes1 * nmodes2, nquad0, nmodes0, 1.0,
236 &inarray[0], nmodes0, base0.data(), nquad0, 0.0, &wsp[0],
237 nmodes1 * nmodes2);
238 Blas::Dgemm('T', 'T', nquad0 * nmodes2, nquad1, nmodes1, 1.0, &wsp[0],
239 nmodes1, base1.data(), nquad1, 0.0, &wsp2[0],
240 nquad0 * nmodes2);
241 Blas::Dgemm('T', 'T', nquad0 * nquad1, nquad2, nmodes2, 1.0, &wsp2[0],
242 nmodes2, base2.data(), nquad2, 0.0, &outarray[0],
243 nquad0 * nquad1);
244 }
245}
static void Dgemm(const char &transa, const char &transb, const int &m, const int &n, const int &k, const double &alpha, const double *a, const int &lda, const double *b, const int &ldb, const double &beta, double *c, const int &ldc)
BLAS level 3: Matrix-matrix multiply C = A x B where op(A)[m x k], op(B)[k x n], C[m x n] DGEMM perfo...
Definition Blas.hpp:383

References ASSERTL1, Blas::Dgemm(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, and Vmath::Vcopy().

◆ v_CalcNumberOfCoefficients()

int Nektar::StdRegions::StdHexExp::v_CalcNumberOfCoefficients ( const std::vector< unsigned int > &  nummodes,
int &  modes_offset 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 715 of file StdHexExp.cpp.

717{
718 int nmodes = nummodes[modes_offset] * nummodes[modes_offset + 1] *
719 nummodes[modes_offset + 2];
720 modes_offset += 3;
721
722 return nmodes;
723}

◆ v_CreateStdMatrix()

DNekMatSharedPtr Nektar::StdRegions::StdHexExp::v_CreateStdMatrix ( const StdMatrixKey mkey)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2258 of file StdHexExp.cpp.

2259{
2260 return v_GenMatrix(mkey);
2261}
DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override

References v_GenMatrix().

◆ v_DetShapeType()

LibUtilities::ShapeType Nektar::StdRegions::StdHexExp::v_DetShapeType ( ) const
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 599 of file StdHexExp.cpp.

References Nektar::LibUtilities::eHexahedron.

◆ v_ExponentialFilter()

void Nektar::StdRegions::StdHexExp::v_ExponentialFilter ( Array< OneD, NekDouble > &  array,
const NekDouble  alpha,
const NekDouble  exponent,
const NekDouble  cutoff 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2461 of file StdHexExp.cpp.

2465{
2466 // Generate an orthogonal expansion
2467 int qa = m_base[0]->GetNumPoints();
2468 int qb = m_base[1]->GetNumPoints();
2469 int qc = m_base[2]->GetNumPoints();
2470 int nmodesA = m_base[0]->GetNumModes();
2471 int nmodesB = m_base[1]->GetNumModes();
2472 int nmodesC = m_base[2]->GetNumModes();
2473 int P = nmodesA - 1;
2474 int Q = nmodesB - 1;
2475 int R = nmodesC - 1;
2476
2477 // Declare orthogonal basis.
2478 LibUtilities::PointsKey pa(qa, m_base[0]->GetPointsType());
2479 LibUtilities::PointsKey pb(qb, m_base[1]->GetPointsType());
2480 LibUtilities::PointsKey pc(qc, m_base[2]->GetPointsType());
2481
2482 LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A, nmodesA, pa);
2483 LibUtilities::BasisKey Bb(LibUtilities::eOrtho_A, nmodesB, pb);
2484 LibUtilities::BasisKey Bc(LibUtilities::eOrtho_A, nmodesC, pc);
2485 StdHexExp OrthoExp(Ba, Bb, Bc);
2486
2487 // Cutoff
2488 int Pcut = cutoff * P;
2489 int Qcut = cutoff * Q;
2490 int Rcut = cutoff * R;
2491
2492 // Project onto orthogonal space.
2493 Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
2494 OrthoExp.FwdTrans(array, orthocoeffs);
2495
2496 //
2497 NekDouble fac, fac1, fac2, fac3;
2498 int index = 0;
2499 for (int i = 0; i < nmodesA; ++i)
2500 {
2501 for (int j = 0; j < nmodesB; ++j)
2502 {
2503 for (int k = 0; k < nmodesC; ++k, ++index)
2504 {
2505 // to filter out only the "high-modes"
2506 if (i > Pcut || j > Qcut || k > Rcut)
2507 {
2508 fac1 = (NekDouble)(i - Pcut) / ((NekDouble)(P - Pcut));
2509 fac2 = (NekDouble)(j - Qcut) / ((NekDouble)(Q - Qcut));
2510 fac3 = (NekDouble)(k - Rcut) / ((NekDouble)(R - Rcut));
2511 fac = max(max(fac1, fac2), fac3);
2512 fac = pow(fac, exponent);
2513 orthocoeffs[index] *= exp(-alpha * fac);
2514 }
2515 }
2516 }
2517 }
2518
2519 // backward transform to physical space
2520 OrthoExp.BwdTrans(orthocoeffs, array);
2521}
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
StdHexExp(const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
Definition StdHexExp.cpp:46
@ P
Monomial polynomials .
Definition BasisType.h:62
@ eOrtho_A
Principle Orthogonal Functions .
Definition BasisType.h:42
scalarT< T > max(scalarT< T > lhs, scalarT< T > rhs)
Definition scalar.hpp:305

References Nektar::StdRegions::StdExpansion::BwdTrans(), Nektar::LibUtilities::eOrtho_A, Nektar::StdRegions::StdExpansion::FwdTrans(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::StdRegions::StdExpansion::m_base, tinysimd::max(), and Nektar::LibUtilities::P.

◆ v_FillMode()

void Nektar::StdRegions::StdHexExp::v_FillMode ( const int  mode,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual
Note
for hexahedral expansions _base[0] (i.e. p) modes run fastest.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 518 of file StdHexExp.cpp.

519{
520 int nquad0 = m_base[0]->GetNumPoints();
521 int nquad1 = m_base[1]->GetNumPoints();
522 int nquad2 = m_base[2]->GetNumPoints();
523
524 Array<OneD, const NekDouble> base0 = m_base[0]->GetBdata();
525 Array<OneD, const NekDouble> base1 = m_base[1]->GetBdata();
526 Array<OneD, const NekDouble> base2 = m_base[2]->GetBdata();
527
528 int btmp0 = m_base[0]->GetNumModes();
529 int btmp1 = m_base[1]->GetNumModes();
530 int mode2 = mode / (btmp0 * btmp1);
531 int mode1 = (mode - mode2 * btmp0 * btmp1) / btmp0;
532 int mode0 = (mode - mode2 * btmp0 * btmp1) % btmp0;
533
534 ASSERTL2(mode == mode2 * btmp0 * btmp1 + mode1 * btmp0 + mode0,
535 "Mode lookup failed.");
536 ASSERTL2(mode < m_ncoeffs,
537 "Calling argument mode is larger than total expansion "
538 "order");
539
540 for (int i = 0; i < nquad1 * nquad2; ++i)
541 {
542 Vmath::Vcopy(nquad0, (NekDouble *)(base0.data() + mode0 * nquad0), 1,
543 &outarray[0] + i * nquad0, 1);
544 }
545
546 for (int j = 0; j < nquad2; ++j)
547 {
548 for (int i = 0; i < nquad0; ++i)
549 {
550 Vmath::Vmul(nquad1, (NekDouble *)(base1.data() + mode1 * nquad1), 1,
551 &outarray[0] + i + j * nquad0 * nquad1, nquad0,
552 &outarray[0] + i + j * nquad0 * nquad1, nquad0);
553 }
554 }
555
556 for (int i = 0; i < nquad2; i++)
557 {
558 Blas::Dscal(nquad0 * nquad1, base2[mode2 * nquad2 + i],
559 &outarray[0] + i * nquad0 * nquad1, 1);
560 }
561}
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition Blas.hpp:149
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition Vmath.hpp:72

References ASSERTL2, Blas::Dscal(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Vmath::Vcopy(), and Vmath::Vmul().

◆ v_FwdTrans()

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

Solves the system \( \mathbf{B^{\top}WB\hat{u}}=\mathbf{B^{\top}Wu^{\delta}} \)

Parameters
inarrayarray of physical quadrature points to be transformed, \( \mathbf{u^{\delta}} \).
outarrayarray of expansion coefficients, \( \mathbf{\hat{u}} \).

Implements Nektar::StdRegions::StdExpansion.

Definition at line 256 of file StdHexExp.cpp.

258{
259 // If using collocation expansion, coefficients match physical
260 // data points so just do a direct copy.
261 if ((m_base[0]->Collocation()) && (m_base[1]->Collocation()) &&
262 (m_base[2]->Collocation()))
263 {
264 Vmath::Vcopy(GetNcoeffs(), &inarray[0], 1, &outarray[0], 1);
265 }
266 else
267 {
268 // Compute B^TWu
269 IProductWRTBase(inarray, outarray);
270
271 // get Mass matrix inverse
272 StdMatrixKey masskey(eInvMass, DetShapeType(), *this);
273 DNekMatSharedPtr matsys = GetStdMatrix(masskey);
274
275 // copy inarray in case inarray == outarray
276 DNekVec in(m_ncoeffs, outarray);
277 DNekVec out(m_ncoeffs, outarray, eWrapper);
278
279 // Solve for coefficients.
280 out = (*matsys) * in;
281 }
282}
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
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...
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
NekVector< NekDouble > DNekVec
std::shared_ptr< DNekMat > DNekMatSharedPtr

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

◆ v_GenMatrix()

DNekMatSharedPtr Nektar::StdRegions::StdHexExp::v_GenMatrix ( const StdMatrixKey mkey)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2170 of file StdHexExp.cpp.

2171{
2172
2173 MatrixType mtype = mkey.GetMatrixType();
2174
2175 DNekMatSharedPtr Mat;
2176
2177 switch (mtype)
2178 {
2180 {
2181 int nq0 = m_base[0]->GetNumPoints();
2182 int nq1 = m_base[1]->GetNumPoints();
2183 int nq2 = m_base[2]->GetNumPoints();
2184 int nq;
2185
2186 // take definition from key
2187 if (mkey.ConstFactorExists(eFactorConst))
2188 {
2189 nq = (int)mkey.GetConstFactor(eFactorConst);
2190 }
2191 else
2192 {
2193 nq = max(nq0, max(nq1, nq2));
2194 }
2195
2196 int neq =
2198 Array<OneD, Array<OneD, NekDouble>> coords(neq);
2199 Array<OneD, NekDouble> coll(3);
2200 Array<OneD, DNekMatSharedPtr> I(3);
2201 Array<OneD, NekDouble> tmp(nq0);
2202
2203 Mat =
2204 MemoryManager<DNekMat>::AllocateSharedPtr(neq, nq0 * nq1 * nq2);
2205 int cnt = 0;
2206
2207 for (int i = 0; i < nq; ++i)
2208 {
2209 for (int j = 0; j < nq; ++j)
2210 {
2211 for (int k = 0; k < nq; ++k, ++cnt)
2212 {
2213 coords[cnt] = Array<OneD, NekDouble>(3);
2214 coords[cnt][0] = -1.0 + 2 * k / (NekDouble)(nq - 1);
2215 coords[cnt][1] = -1.0 + 2 * j / (NekDouble)(nq - 1);
2216 coords[cnt][2] = -1.0 + 2 * i / (NekDouble)(nq - 1);
2217 }
2218 }
2219 }
2220
2221 for (int i = 0; i < neq; ++i)
2222 {
2223 LocCoordToLocCollapsed(coords[i], coll);
2224
2225 I[0] = m_base[0]->GetI(coll);
2226 I[1] = m_base[1]->GetI(coll + 1);
2227 I[2] = m_base[2]->GetI(coll + 2);
2228
2229 // interpolate first coordinate direction
2230 NekDouble fac;
2231 for (int k = 0; k < nq2; ++k)
2232 {
2233 for (int j = 0; j < nq1; ++j)
2234 {
2235
2236 fac = (I[1]->GetPtr())[j] * (I[2]->GetPtr())[k];
2237 Vmath::Smul(nq0, fac, I[0]->GetPtr(), 1, tmp, 1);
2238
2239 Vmath::Vcopy(nq0, &tmp[0], 1,
2240 Mat->GetRawPtr() + k * nq0 * nq1 * neq +
2241 j * nq0 * neq + i,
2242 neq);
2243 }
2244 }
2245 }
2246 }
2247 break;
2248 default:
2249 {
2251 }
2252 break;
2253 }
2254
2255 return Mat;
2256}
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
int getNumberOfCoefficients(int Na, int Nb, int Nc)
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

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::StdRegions::StdMatrixKey::ConstFactorExists(), Nektar::StdRegions::StdExpansion::CreateGeneralMatrix(), Nektar::StdRegions::eFactorConst, Nektar::StdRegions::ePhysInterpToEquiSpaced, Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdMatrixKey::GetMatrixType(), Nektar::LibUtilities::StdHexData::getNumberOfCoefficients(), Nektar::StdRegions::StdExpansion::LocCoordToLocCollapsed(), Nektar::StdRegions::StdExpansion::m_base, tinysimd::max(), Vmath::Smul(), and Vmath::Vcopy().

Referenced by v_CreateStdMatrix().

◆ v_GetBoundaryMap()

void Nektar::StdRegions::StdHexExp::v_GetBoundaryMap ( Array< OneD, unsigned int > &  outarray)
overrideprotectedvirtual
Parameters
outarrayStorage for computed map.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1016 of file StdHexExp.cpp.

1017{
1020 "BasisType is not a boundary interior form");
1023 "BasisType is not a boundary interior form");
1026 "BasisType is not a boundary interior form");
1027
1028 int i;
1029 int nummodes[3] = {m_base[0]->GetNumModes(), m_base[1]->GetNumModes(),
1030 m_base[2]->GetNumModes()};
1031
1032 int nBndCoeffs = NumBndryCoeffs();
1033
1034 if (outarray.size() != nBndCoeffs)
1035 {
1036 outarray = Array<OneD, unsigned int>(nBndCoeffs);
1037 }
1038
1039 const LibUtilities::BasisType Btype[3] = {GetBasisType(0), GetBasisType(1),
1040 GetBasisType(2)};
1041
1042 int p, q, r;
1043 int cnt = 0;
1044
1045 int BndIdx[3][2];
1046 int IntIdx[3][2];
1047
1048 for (i = 0; i < 3; i++)
1049 {
1050 BndIdx[i][0] = 0;
1051
1052 if (Btype[i] == LibUtilities::eModified_A)
1053 {
1054 BndIdx[i][1] = 1;
1055 IntIdx[i][0] = 2;
1056 IntIdx[i][1] = nummodes[i];
1057 }
1058 else
1059 {
1060 BndIdx[i][1] = nummodes[i] - 1;
1061 IntIdx[i][0] = 1;
1062 IntIdx[i][1] = nummodes[i] - 1;
1063 }
1064 }
1065
1066 for (i = 0; i < 2; i++)
1067 {
1068 r = BndIdx[2][i];
1069 for (q = 0; q < nummodes[1]; q++)
1070 {
1071 for (p = 0; p < nummodes[0]; p++)
1072 {
1073 outarray[cnt++] =
1074 r * nummodes[0] * nummodes[1] + q * nummodes[0] + p;
1075 }
1076 }
1077 }
1078
1079 for (r = IntIdx[2][0]; r < IntIdx[2][1]; r++)
1080 {
1081 for (i = 0; i < 2; i++)
1082 {
1083 q = BndIdx[1][i];
1084 for (p = 0; p < nummodes[0]; p++)
1085 {
1086 outarray[cnt++] =
1087 r * nummodes[0] * nummodes[1] + q * nummodes[0] + p;
1088 }
1089 }
1090
1091 for (q = IntIdx[1][0]; q < IntIdx[1][1]; q++)
1092 {
1093 for (i = 0; i < 2; i++)
1094 {
1095 p = BndIdx[0][i];
1096 outarray[cnt++] =
1097 r * nummodes[0] * nummodes[1] + q * nummodes[0] + p;
1098 }
1099 }
1100 }
1101
1102 sort(outarray.data(), outarray.data() + nBndCoeffs);
1103}
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition BasisType.h:56
@ eModified_A
Principle Modified Functions .
Definition BasisType.h:48
std::vector< double > p(NPUPPER)
std::vector< double > q(NPUPPER *NPUPPER)

References ASSERTL1, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::m_base, and Nektar::StdRegions::StdExpansion::NumBndryCoeffs().

◆ v_GetCoords()

void Nektar::StdRegions::StdHexExp::v_GetCoords ( Array< OneD, NekDouble > &  coords_x,
Array< OneD, NekDouble > &  coords_y,
Array< OneD, NekDouble > &  coords_z 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 751 of file StdHexExp.cpp.

754{
755 Array<OneD, const NekDouble> eta_x = m_base[0]->GetZ();
756 Array<OneD, const NekDouble> eta_y = m_base[1]->GetZ();
757 Array<OneD, const NekDouble> eta_z = m_base[2]->GetZ();
758 int Qx = GetNumPoints(0);
759 int Qy = GetNumPoints(1);
760 int Qz = GetNumPoints(2);
761
762 // Convert collapsed coordinates into cartesian coordinates:
763 // eta --> xi
764 for (int k = 0; k < Qz; ++k)
765 {
766 for (int j = 0; j < Qy; ++j)
767 {
768 for (int i = 0; i < Qx; ++i)
769 {
770 int s = i + Qx * (j + Qy * k);
771 xi_x[s] = eta_x[i];
772 xi_y[s] = eta_y[j];
773 xi_z[s] = eta_z[k];
774 }
775 }
776 }
777}

References Nektar::StdRegions::StdExpansion::GetNumPoints(), and Nektar::StdRegions::StdExpansion::m_base.

◆ v_GetEdgeInteriorToElementMap()

void Nektar::StdRegions::StdHexExp::v_GetEdgeInteriorToElementMap ( const int  eid,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray,
const Orientation  edgeOrient = eDir1FwdDir1_Dir2FwdDir2 
)
overrideprotectedvirtual
Parameters
eidThe edge to compute the numbering for.
edgeOrientOrientation of the edge.
maparrayStorage for computed mapping array.
signarray?

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 1535 of file StdHexExp.cpp.

1538{
1541 "BasisType is not a boundary interior form");
1544 "BasisType is not a boundary interior form");
1547 "BasisType is not a boundary interior form");
1548
1549 ASSERTL1((eid >= 0) && (eid < 12),
1550 "local edge id must be between 0 and 11");
1551
1552 int nEdgeIntCoeffs = GetEdgeNcoeffs(eid) - 2;
1553
1554 if (maparray.size() != nEdgeIntCoeffs)
1555 {
1556 maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
1557 }
1558
1559 if (signarray.size() != nEdgeIntCoeffs)
1560 {
1561 signarray = Array<OneD, int>(nEdgeIntCoeffs, 1);
1562 }
1563 else
1564 {
1565 fill(signarray.data(), signarray.data() + nEdgeIntCoeffs, 1);
1566 }
1567
1568 int nummodes[3] = {m_base[0]->GetNumModes(), m_base[1]->GetNumModes(),
1569 m_base[2]->GetNumModes()};
1570
1571 const LibUtilities::BasisType bType[3] = {GetBasisType(0), GetBasisType(1),
1572 GetBasisType(2)};
1573
1574 bool reverseOrdering = false;
1575 bool signChange = false;
1576
1577 int IdxRange[3][2] = {{0, 0}, {0, 0}, {0, 0}};
1578
1579 switch (eid)
1580 {
1581 case 0:
1582 case 1:
1583 case 2:
1584 case 3:
1585 {
1586 IdxRange[2][0] = 0;
1587 IdxRange[2][1] = 1;
1588 }
1589 break;
1590 case 8:
1591 case 9:
1592 case 10:
1593 case 11:
1594 {
1595 if (bType[2] == LibUtilities::eGLL_Lagrange)
1596 {
1597 IdxRange[2][0] = nummodes[2] - 1;
1598 IdxRange[2][1] = nummodes[2];
1599 }
1600 else
1601 {
1602 IdxRange[2][0] = 1;
1603 IdxRange[2][1] = 2;
1604 }
1605 }
1606 break;
1607 case 4:
1608 case 5:
1609 case 6:
1610 case 7:
1611 {
1612 if (bType[2] == LibUtilities::eGLL_Lagrange)
1613 {
1614 IdxRange[2][0] = 1;
1615 IdxRange[2][1] = nummodes[2] - 1;
1616
1617 if (edgeOrient == eBackwards)
1618 {
1619 reverseOrdering = true;
1620 }
1621 }
1622 else
1623 {
1624 IdxRange[2][0] = 2;
1625 IdxRange[2][1] = nummodes[2];
1626
1627 if (edgeOrient == eBackwards)
1628 {
1629 signChange = true;
1630 }
1631 }
1632 }
1633 break;
1634 }
1635
1636 switch (eid)
1637 {
1638 case 0:
1639 case 4:
1640 case 5:
1641 case 8:
1642 {
1643 IdxRange[1][0] = 0;
1644 IdxRange[1][1] = 1;
1645 }
1646 break;
1647 case 2:
1648 case 6:
1649 case 7:
1650 case 10:
1651 {
1652 if (bType[1] == LibUtilities::eGLL_Lagrange)
1653 {
1654 IdxRange[1][0] = nummodes[1] - 1;
1655 IdxRange[1][1] = nummodes[1];
1656 }
1657 else
1658 {
1659 IdxRange[1][0] = 1;
1660 IdxRange[1][1] = 2;
1661 }
1662 }
1663 break;
1664 case 1:
1665 case 9:
1666 {
1667 if (bType[1] == LibUtilities::eGLL_Lagrange)
1668 {
1669 IdxRange[1][0] = 1;
1670 IdxRange[1][1] = nummodes[1] - 1;
1671
1672 if (edgeOrient == eBackwards)
1673 {
1674 reverseOrdering = true;
1675 }
1676 }
1677 else
1678 {
1679 IdxRange[1][0] = 2;
1680 IdxRange[1][1] = nummodes[1];
1681
1682 if (edgeOrient == eBackwards)
1683 {
1684 signChange = true;
1685 }
1686 }
1687 }
1688 break;
1689 case 3:
1690 case 11:
1691 {
1692 if (bType[1] == LibUtilities::eGLL_Lagrange)
1693 {
1694 IdxRange[1][0] = 1;
1695 IdxRange[1][1] = nummodes[1] - 1;
1696
1697 if (edgeOrient == eForwards)
1698 {
1699 reverseOrdering = true;
1700 }
1701 }
1702 else
1703 {
1704 IdxRange[1][0] = 2;
1705 IdxRange[1][1] = nummodes[1];
1706
1707 if (edgeOrient == eForwards)
1708 {
1709 signChange = true;
1710 }
1711 }
1712 }
1713 break;
1714 }
1715
1716 switch (eid)
1717 {
1718 case 3:
1719 case 4:
1720 case 7:
1721 case 11:
1722 {
1723 IdxRange[0][0] = 0;
1724 IdxRange[0][1] = 1;
1725 }
1726 break;
1727 case 1:
1728 case 5:
1729 case 6:
1730 case 9:
1731 {
1732 if (bType[0] == LibUtilities::eGLL_Lagrange)
1733 {
1734 IdxRange[0][0] = nummodes[0] - 1;
1735 IdxRange[0][1] = nummodes[0];
1736 }
1737 else
1738 {
1739 IdxRange[0][0] = 1;
1740 IdxRange[0][1] = 2;
1741 }
1742 }
1743 break;
1744 case 0:
1745 case 8:
1746 {
1747 if (bType[0] == LibUtilities::eGLL_Lagrange)
1748 {
1749 IdxRange[0][0] = 1;
1750 IdxRange[0][1] = nummodes[0] - 1;
1751
1752 if (edgeOrient == eBackwards)
1753 {
1754 reverseOrdering = true;
1755 }
1756 }
1757 else
1758 {
1759 IdxRange[0][0] = 2;
1760 IdxRange[0][1] = nummodes[0];
1761
1762 if (edgeOrient == eBackwards)
1763 {
1764 signChange = true;
1765 }
1766 }
1767 }
1768 break;
1769 case 2:
1770 case 10:
1771 {
1772 if (bType[0] == LibUtilities::eGLL_Lagrange)
1773 {
1774 IdxRange[0][0] = 1;
1775 IdxRange[0][1] = nummodes[0] - 1;
1776
1777 if (edgeOrient == eForwards)
1778 {
1779 reverseOrdering = true;
1780 }
1781 }
1782 else
1783 {
1784 IdxRange[0][0] = 2;
1785 IdxRange[0][1] = nummodes[0];
1786
1787 if (edgeOrient == eForwards)
1788 {
1789 signChange = true;
1790 }
1791 }
1792 }
1793 break;
1794 }
1795
1796 int cnt = 0;
1797
1798 for (int r = IdxRange[2][0]; r < IdxRange[2][1]; r++)
1799 {
1800 for (int q = IdxRange[1][0]; q < IdxRange[1][1]; q++)
1801 {
1802 for (int p = IdxRange[0][0]; p < IdxRange[0][1]; p++)
1803 {
1804 maparray[cnt++] =
1805 r * nummodes[0] * nummodes[1] + q * nummodes[0] + p;
1806 }
1807 }
1808 }
1809
1810 if (reverseOrdering)
1811 {
1812 reverse(maparray.data(), maparray.data() + nEdgeIntCoeffs);
1813 }
1814
1815 if (signChange)
1816 {
1817 for (int p = 1; p < nEdgeIntCoeffs; p += 2)
1818 {
1819 signarray[p] = -1;
1820 }
1821 }
1822}
int GetEdgeNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th edge.

References ASSERTL1, Nektar::StdRegions::eBackwards, Nektar::StdRegions::eForwards, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion3D::GetEdgeNcoeffs(), and Nektar::StdRegions::StdExpansion::m_base.

◆ v_GetEdgeNcoeffs()

int Nektar::StdRegions::StdHexExp::v_GetEdgeNcoeffs ( const int  i) const
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 2152 of file StdHexExp.cpp.

2153{
2154 ASSERTL2((i >= 0) && (i <= 11), "edge id is out of range");
2155
2156 if ((i == 0) || (i == 2) || (i == 8) || (i == 10))
2157 {
2158 return GetBasisNumModes(0);
2159 }
2160 else if ((i == 1) || (i == 3) || (i == 9) || (i == 11))
2161 {
2162 return GetBasisNumModes(1);
2163 }
2164 else
2165 {
2166 return GetBasisNumModes(2);
2167 }
2168}
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.

References ASSERTL2, and Nektar::StdRegions::StdExpansion::GetBasisNumModes().

◆ v_GetElmtTraceToTraceMap()

void Nektar::StdRegions::StdHexExp::v_GetElmtTraceToTraceMap ( const unsigned int  fid,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray,
Orientation  faceOrient,
int  P,
int  Q 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1240 of file StdHexExp.cpp.

1244{
1245 int i, j;
1246 int nummodesA = 0, nummodesB = 0;
1247
1249 GetBasisType(0) == GetBasisType(2),
1250 "Method only implemented if BasisType is indentical in "
1251 "all directions");
1254 "Method only implemented for Modified_A or "
1255 "GLL_Lagrange BasisType");
1256
1257 const int nummodes0 = m_base[0]->GetNumModes();
1258 const int nummodes1 = m_base[1]->GetNumModes();
1259 const int nummodes2 = m_base[2]->GetNumModes();
1260
1261 switch (fid)
1262 {
1263 case 0:
1264 case 5:
1265 nummodesA = nummodes0;
1266 nummodesB = nummodes1;
1267 break;
1268 case 1:
1269 case 3:
1270 nummodesA = nummodes0;
1271 nummodesB = nummodes2;
1272 break;
1273 case 2:
1274 case 4:
1275 nummodesA = nummodes1;
1276 nummodesB = nummodes2;
1277 break;
1278 default:
1279 ASSERTL0(false, "fid must be between 0 and 5");
1280 }
1281
1282 if (P == -1)
1283 {
1284 P = nummodesA;
1285 Q = nummodesB;
1286 }
1287
1288 bool modified = (GetBasisType(0) == LibUtilities::eModified_A);
1289
1290 // check that
1291 if (modified == false)
1292 {
1293 ASSERTL1((P == nummodesA) || (Q == nummodesB),
1294 "Different trace space face dimention "
1295 "and element face dimention not possible for "
1296 "GLL-Lagrange bases");
1297 }
1298
1299 int nFaceCoeffs = P * Q;
1300
1301 if (maparray.size() != nFaceCoeffs)
1302 {
1303 maparray = Array<OneD, unsigned int>(nFaceCoeffs);
1304 }
1305
1306 // fill default mapping as increasing index
1307 for (i = 0; i < nFaceCoeffs; ++i)
1308 {
1309 maparray[i] = i;
1310 }
1311
1312 if (signarray.size() != nFaceCoeffs)
1313 {
1314 signarray = Array<OneD, int>(nFaceCoeffs, 1);
1315 }
1316 else
1317 {
1318 fill(signarray.data(), signarray.data() + nFaceCoeffs, 1);
1319 }
1320
1321 // setup indexing to manage transpose directions
1322 Array<OneD, int> arrayindx(nFaceCoeffs);
1323 for (i = 0; i < Q; i++)
1324 {
1325 for (j = 0; j < P; j++)
1326 {
1327 if (faceOrient < eDir1FwdDir2_Dir2FwdDir1)
1328 {
1329 arrayindx[i * P + j] = i * P + j;
1330 }
1331 else
1332 {
1333 arrayindx[i * P + j] = j * Q + i;
1334 }
1335 }
1336 }
1337
1338 // zero signmap and set maparray to zero if elemental
1339 // modes are not as large as face models
1340 for (i = 0; i < nummodesB; i++)
1341 {
1342 for (j = nummodesA; j < P; j++)
1343 {
1344 signarray[arrayindx[i * P + j]] = 0.0;
1345 maparray[arrayindx[i * P + j]] = maparray[0];
1346 }
1347 }
1348
1349 for (i = nummodesB; i < Q; i++)
1350 {
1351 for (j = 0; j < P; j++)
1352 {
1353 signarray[arrayindx[i * P + j]] = 0.0;
1354 maparray[arrayindx[i * P + j]] = maparray[0];
1355 }
1356 }
1357
1358 // zero signmap and set maparray to zero entry if
1359 // elemental modes are not as large as face modesl
1360 for (i = 0; i < Q; i++)
1361 {
1362 // fill values into map array of trace size
1363 // for element face index
1364 for (j = 0; j < P; j++)
1365 {
1366 maparray[arrayindx[i * P + j]] = i * nummodesA + j;
1367 }
1368
1369 // zero values if P > numModesA
1370 for (j = nummodesA; j < P; j++)
1371 {
1372 signarray[arrayindx[i * P + j]] = 0.0;
1373 maparray[arrayindx[i * P + j]] = maparray[0];
1374 }
1375 }
1376
1377 // zero values if Q > numModesB
1378 for (i = nummodesB; i < Q; i++)
1379 {
1380 for (j = 0; j < P; j++)
1381 {
1382 signarray[arrayindx[i * P + j]] = 0.0;
1383 maparray[arrayindx[i * P + j]] = maparray[0];
1384 }
1385 }
1386
1387 // Now reorientate indices accordign to orientation
1388 if ((faceOrient == eDir1FwdDir1_Dir2BwdDir2) ||
1389 (faceOrient == eDir1BwdDir1_Dir2BwdDir2) ||
1390 (faceOrient == eDir1BwdDir2_Dir2FwdDir1) ||
1391 (faceOrient == eDir1BwdDir2_Dir2BwdDir1))
1392 {
1393 if (faceOrient < eDir1FwdDir2_Dir2FwdDir1)
1394 {
1395 if (modified)
1396 {
1397 for (i = 3; i < Q; i += 2)
1398 {
1399 for (j = 0; j < P; j++)
1400 {
1401 signarray[arrayindx[i * P + j]] *= -1;
1402 }
1403 }
1404
1405 for (i = 0; i < P; i++)
1406 {
1407 swap(maparray[i], maparray[i + P]);
1408 swap(signarray[i], signarray[i + P]);
1409 }
1410 }
1411 else
1412 {
1413 for (i = 0; i < P; i++)
1414 {
1415 for (j = 0; j < Q / 2; j++)
1416 {
1417 swap(maparray[i + j * P],
1418 maparray[i + P * Q - P - j * P]);
1419 swap(signarray[i + j * P],
1420 signarray[i + P * Q - P - j * P]);
1421 }
1422 }
1423 }
1424 }
1425 else
1426 {
1427 if (modified)
1428 {
1429 for (i = 0; i < Q; i++)
1430 {
1431 for (j = 3; j < P; j += 2)
1432 {
1433 signarray[arrayindx[i * P + j]] *= -1;
1434 }
1435 }
1436
1437 for (i = 0; i < Q; i++)
1438 {
1439 swap(maparray[i], maparray[i + Q]);
1440 swap(signarray[i], signarray[i + Q]);
1441 }
1442 }
1443 else
1444 {
1445 for (i = 0; i < P; i++)
1446 {
1447 for (j = 0; j < Q / 2; j++)
1448 {
1449 swap(maparray[i * Q + j], maparray[i * Q + Q - 1 - j]);
1450 swap(signarray[i * Q + j],
1451 signarray[i * Q + Q - 1 - j]);
1452 }
1453 }
1454 }
1455 }
1456 }
1457
1458 if ((faceOrient == eDir1BwdDir1_Dir2FwdDir2) ||
1459 (faceOrient == eDir1BwdDir1_Dir2BwdDir2) ||
1460 (faceOrient == eDir1FwdDir2_Dir2BwdDir1) ||
1461 (faceOrient == eDir1BwdDir2_Dir2BwdDir1))
1462 {
1463 if (faceOrient < eDir1FwdDir2_Dir2FwdDir1)
1464 {
1465 if (modified)
1466 {
1467 for (i = 0; i < Q; i++)
1468 {
1469 for (j = 3; j < P; j += 2)
1470 {
1471 signarray[arrayindx[i * P + j]] *= -1;
1472 }
1473 }
1474
1475 for (i = 0; i < Q; i++)
1476 {
1477 swap(maparray[i * P], maparray[i * P + 1]);
1478 swap(signarray[i * P], signarray[i * P + 1]);
1479 }
1480 }
1481 else
1482 {
1483 for (i = 0; i < Q; i++)
1484 {
1485 for (j = 0; j < P / 2; j++)
1486 {
1487 swap(maparray[i * P + j], maparray[i * P + P - 1 - j]);
1488 swap(signarray[i * P + j],
1489 signarray[i * P + P - 1 - j]);
1490 }
1491 }
1492 }
1493 }
1494 else
1495 {
1496 if (modified)
1497 {
1498 for (i = 3; i < Q; i += 2)
1499 {
1500 for (j = 0; j < P; j++)
1501 {
1502 signarray[arrayindx[i * P + j]] *= -1;
1503 }
1504 }
1505
1506 for (i = 0; i < P; i++)
1507 {
1508 swap(maparray[i * Q], maparray[i * Q + 1]);
1509 swap(signarray[i * Q], signarray[i * Q + 1]);
1510 }
1511 }
1512 else
1513 {
1514 for (i = 0; i < Q; i++)
1515 {
1516 for (j = 0; j < P / 2; j++)
1517 {
1518 swap(maparray[i + j * Q],
1519 maparray[i + P * Q - Q - j * Q]);
1520 swap(signarray[i + j * Q],
1521 signarray[i + P * Q - Q - j * Q]);
1522 }
1523 }
1524 }
1525 }
1526 }
1527}
#define ASSERTL0(condition, msg)

References ASSERTL0, ASSERTL1, Nektar::StdRegions::eDir1BwdDir1_Dir2BwdDir2, Nektar::StdRegions::eDir1BwdDir1_Dir2FwdDir2, Nektar::StdRegions::eDir1BwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1BwdDir2_Dir2FwdDir1, Nektar::StdRegions::eDir1FwdDir1_Dir2BwdDir2, Nektar::StdRegions::eDir1FwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1FwdDir2_Dir2FwdDir1, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::m_base, and Nektar::LibUtilities::P.

◆ v_GetInteriorMap()

void Nektar::StdRegions::StdHexExp::v_GetInteriorMap ( Array< OneD, unsigned int > &  outarray)
overrideprotectedvirtual
Parameters
outarrayStorage area for computed map.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 955 of file StdHexExp.cpp.

956{
959 "BasisType is not a boundary interior form");
962 "BasisType is not a boundary interior form");
965 "BasisType is not a boundary interior form");
966
967 int i;
968 int nummodes[3] = {m_base[0]->GetNumModes(), m_base[1]->GetNumModes(),
969 m_base[2]->GetNumModes()};
970
971 int nIntCoeffs = m_ncoeffs - NumBndryCoeffs();
972
973 if (outarray.size() != nIntCoeffs)
974 {
975 outarray = Array<OneD, unsigned int>(nIntCoeffs);
976 }
977
978 const LibUtilities::BasisType Btype[3] = {GetBasisType(0), GetBasisType(1),
979 GetBasisType(2)};
980
981 int p, q, r;
982 int cnt = 0;
983
984 int IntIdx[3][2];
985
986 for (i = 0; i < 3; i++)
987 {
988 if (Btype[i] == LibUtilities::eModified_A)
989 {
990 IntIdx[i][0] = 2;
991 IntIdx[i][1] = nummodes[i];
992 }
993 else
994 {
995 IntIdx[i][0] = 1;
996 IntIdx[i][1] = nummodes[i] - 1;
997 }
998 }
999
1000 for (r = IntIdx[2][0]; r < IntIdx[2][1]; r++)
1001 {
1002 for (q = IntIdx[1][0]; q < IntIdx[1][1]; q++)
1003 {
1004 for (p = IntIdx[0][0]; p < IntIdx[0][1]; p++)
1005 {
1006 outarray[cnt++] =
1007 r * nummodes[0] * nummodes[1] + q * nummodes[0] + p;
1008 }
1009 }
1010 }
1011}

References ASSERTL1, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, and Nektar::StdRegions::StdExpansion::NumBndryCoeffs().

◆ v_GetNedges()

int Nektar::StdRegions::StdHexExp::v_GetNedges ( void  ) const
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 589 of file StdHexExp.cpp.

590{
591 return 12;
592}

◆ v_GetNtraces()

int Nektar::StdRegions::StdHexExp::v_GetNtraces ( ) const
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 594 of file StdHexExp.cpp.

595{
596 return 6;
597}

◆ v_GetNverts()

int Nektar::StdRegions::StdHexExp::v_GetNverts ( ) const
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 584 of file StdHexExp.cpp.

585{
586 return 8;
587}

◆ v_GetSimplexEquiSpacedConnectivity()

void Nektar::StdRegions::StdHexExp::v_GetSimplexEquiSpacedConnectivity ( Array< OneD, int > &  conn,
bool  standard = true 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2523 of file StdHexExp.cpp.

2525{
2526 int np0 = m_base[0]->GetNumPoints();
2527 int np1 = m_base[1]->GetNumPoints();
2528 int np2 = m_base[2]->GetNumPoints();
2529 int np = max(np0, max(np1, np2));
2530
2531 conn = Array<OneD, int>(6 * (np - 1) * (np - 1) * (np - 1));
2532
2533 int row = 0;
2534 int rowp1 = 0;
2535 int cnt = 0;
2536 int plane = 0;
2537 for (int i = 0; i < np - 1; ++i)
2538 {
2539 for (int j = 0; j < np - 1; ++j)
2540 {
2541 rowp1 += np;
2542 for (int k = 0; k < np - 1; ++k)
2543 {
2544 conn[cnt++] = plane + row + k;
2545 conn[cnt++] = plane + row + k + 1;
2546 conn[cnt++] = plane + rowp1 + k;
2547
2548 conn[cnt++] = plane + rowp1 + k + 1;
2549 conn[cnt++] = plane + rowp1 + k;
2550 conn[cnt++] = plane + row + k + 1;
2551 }
2552 row += np;
2553 }
2554 plane += np * np;
2555 }
2556}

References Nektar::StdRegions::StdExpansion::m_base, and tinysimd::max().

◆ v_GetTraceBasisKey()

const LibUtilities::BasisKey Nektar::StdRegions::StdHexExp::v_GetTraceBasisKey ( const int  i,
const int  k,
bool  useGLL = false 
) const
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 725 of file StdHexExp.cpp.

727{
728 ASSERTL2(i >= 0 && i <= 5, "face id is out of range");
729 ASSERTL2(k >= 0 && k <= 1, "basis key id is out of range");
730
731 int dir = k;
732 switch (i)
733 {
734 case 0:
735 case 5:
736 dir = k;
737 break;
738 case 1:
739 case 3:
740 dir = 2 * k;
741 break;
742 case 2:
743 case 4:
744 dir = k + 1;
745 break;
746 }
747
748 return EvaluateQuadFaceBasisKey(k, m_base[dir]);
749}
LibUtilities::BasisKey EvaluateQuadFaceBasisKey(const int facedir, const LibUtilities::BasisSharedPtr &faceDirBasis)

References ASSERTL2, Nektar::StdRegions::EvaluateQuadFaceBasisKey(), and Nektar::StdRegions::StdExpansion::m_base.

◆ v_GetTraceCoeffMap()

void Nektar::StdRegions::StdHexExp::v_GetTraceCoeffMap ( const unsigned int  fid,
Array< OneD, unsigned int > &  maparray 
)
overrideprotectedvirtual

Only for basis type Modified_A or GLL_LAGRANGE in all directions.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1116 of file StdHexExp.cpp.

1118{
1119 int i, j;
1120 int nummodesA = 0, nummodesB = 0;
1121
1123 GetBasisType(0) == GetBasisType(2),
1124 "Method only implemented if BasisType is indentical in "
1125 "all directions");
1128 "Method only implemented for Modified_A or "
1129 "GLL_Lagrange BasisType");
1130
1131 const int nummodes0 = m_base[0]->GetNumModes();
1132 const int nummodes1 = m_base[1]->GetNumModes();
1133 const int nummodes2 = m_base[2]->GetNumModes();
1134
1135 switch (fid)
1136 {
1137 case 0:
1138 case 5:
1139 nummodesA = nummodes0;
1140 nummodesB = nummodes1;
1141 break;
1142 case 1:
1143 case 3:
1144 nummodesA = nummodes0;
1145 nummodesB = nummodes2;
1146 break;
1147 case 2:
1148 case 4:
1149 nummodesA = nummodes1;
1150 nummodesB = nummodes2;
1151 break;
1152 default:
1153 ASSERTL0(false, "fid must be between 0 and 5");
1154 }
1155
1156 int nFaceCoeffs = nummodesA * nummodesB;
1157
1158 if (maparray.size() != nFaceCoeffs)
1159 {
1160 maparray = Array<OneD, unsigned int>(nFaceCoeffs);
1161 }
1162
1163 bool modified = (GetBasisType(0) == LibUtilities::eModified_A);
1164
1165 int offset = 0;
1166 int jump1 = 1;
1167 int jump2 = 1;
1168
1169 switch (fid)
1170 {
1171 case 5:
1172 {
1173 if (modified)
1174 {
1175 offset = nummodes0 * nummodes1;
1176 }
1177 else
1178 {
1179 offset = (nummodes2 - 1) * nummodes0 * nummodes1;
1180 jump1 = nummodes0;
1181 }
1182 }
1183 /* Falls through. */
1184 case 0:
1185 {
1186 jump1 = nummodes0;
1187 break;
1188 }
1189 case 3:
1190 {
1191 if (modified)
1192 {
1193 offset = nummodes0;
1194 }
1195 else
1196 {
1197 offset = nummodes0 * (nummodes1 - 1);
1198 jump1 = nummodes0 * nummodes1;
1199 }
1200 }
1201 /* Falls through. */
1202 case 1:
1203 {
1204 jump1 = nummodes0 * nummodes1;
1205 break;
1206 }
1207 case 2:
1208 {
1209 if (modified)
1210 {
1211 offset = 1;
1212 }
1213 else
1214 {
1215 offset = nummodes0 - 1;
1216 jump1 = nummodes0 * nummodes1;
1217 jump2 = nummodes0;
1218 }
1219 }
1220 /* Falls through. */
1221 case 4:
1222 {
1223 jump1 = nummodes0 * nummodes1;
1224 jump2 = nummodes0;
1225 break;
1226 }
1227 default:
1228 ASSERTL0(false, "fid must be between 0 and 5");
1229 }
1230
1231 for (i = 0; i < nummodesB; i++)
1232 {
1233 for (j = 0; j < nummodesA; j++)
1234 {
1235 maparray[i * nummodesA + j] = i * jump1 + j * jump2 + offset;
1236 }
1237 }
1238}

References ASSERTL0, ASSERTL1, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), and Nektar::StdRegions::StdExpansion::m_base.

◆ v_GetTraceInteriorToElementMap()

void Nektar::StdRegions::StdHexExp::v_GetTraceInteriorToElementMap ( const int  fid,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray,
const Orientation  faceOrient = eDir1FwdDir1_Dir2FwdDir2 
)
overrideprotectedvirtual

Generate mapping describing which elemental modes lie on the interior of a given face. Accounts for face orientation.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1828 of file StdHexExp.cpp.

1831{
1834 "BasisType is not a boundary interior form");
1837 "BasisType is not a boundary interior form");
1840 "BasisType is not a boundary interior form");
1841
1842 ASSERTL1((fid >= 0) && (fid < 6), "local face id must be between 0 and 5");
1843
1844 int nFaceIntCoeffs = v_GetTraceIntNcoeffs(fid);
1845
1846 if (maparray.size() != nFaceIntCoeffs)
1847 {
1848 maparray = Array<OneD, unsigned int>(nFaceIntCoeffs);
1849 }
1850
1851 if (signarray.size() != nFaceIntCoeffs)
1852 {
1853 signarray = Array<OneD, int>(nFaceIntCoeffs, 1);
1854 }
1855 else
1856 {
1857 fill(signarray.data(), signarray.data() + nFaceIntCoeffs, 1);
1858 }
1859
1860 int nummodes[3] = {m_base[0]->GetNumModes(), m_base[1]->GetNumModes(),
1861 m_base[2]->GetNumModes()};
1862
1863 const LibUtilities::BasisType bType[3] = {GetBasisType(0), GetBasisType(1),
1864 GetBasisType(2)};
1865
1866 int nummodesA = 0;
1867 int nummodesB = 0;
1868
1869 // Determine the number of modes in face directions A & B based
1870 // on the face index given.
1871 switch (fid)
1872 {
1873 case 0:
1874 case 5:
1875 {
1876 nummodesA = nummodes[0];
1877 nummodesB = nummodes[1];
1878 }
1879 break;
1880 case 1:
1881 case 3:
1882 {
1883 nummodesA = nummodes[0];
1884 nummodesB = nummodes[2];
1885 }
1886 break;
1887 case 2:
1888 case 4:
1889 {
1890 nummodesA = nummodes[1];
1891 nummodesB = nummodes[2];
1892 }
1893 }
1894
1895 Array<OneD, int> arrayindx(nFaceIntCoeffs);
1896
1897 // Create a mapping array to account for transposition of the
1898 // coordinates due to face orientation.
1899 for (int i = 0; i < (nummodesB - 2); i++)
1900 {
1901 for (int j = 0; j < (nummodesA - 2); j++)
1902 {
1903 if (faceOrient < eDir1FwdDir2_Dir2FwdDir1)
1904 {
1905 arrayindx[i * (nummodesA - 2) + j] = i * (nummodesA - 2) + j;
1906 }
1907 else
1908 {
1909 arrayindx[i * (nummodesA - 2) + j] = j * (nummodesB - 2) + i;
1910 }
1911 }
1912 }
1913
1914 int IdxRange[3][2];
1915 int Incr[3];
1916
1917 Array<OneD, int> sign0(nummodes[0], 1);
1918 Array<OneD, int> sign1(nummodes[1], 1);
1919 Array<OneD, int> sign2(nummodes[2], 1);
1920
1921 // Set the upper and lower bounds, and increment for the faces
1922 // involving the first coordinate direction.
1923 switch (fid)
1924 {
1925 case 0: // bottom face
1926 {
1927 IdxRange[2][0] = 0;
1928 IdxRange[2][1] = 1;
1929 Incr[2] = 1;
1930 }
1931 break;
1932 case 5: // top face
1933 {
1934 if (bType[2] == LibUtilities::eGLL_Lagrange)
1935 {
1936 IdxRange[2][0] = nummodes[2] - 1;
1937 IdxRange[2][1] = nummodes[2];
1938 Incr[2] = 1;
1939 }
1940 else
1941 {
1942 IdxRange[2][0] = 1;
1943 IdxRange[2][1] = 2;
1944 Incr[2] = 1;
1945 }
1946 }
1947 break;
1948 default: // all other faces
1949 {
1950 if (bType[2] == LibUtilities::eGLL_Lagrange)
1951 {
1952 if (((int)(faceOrient - eDir1FwdDir1_Dir2FwdDir2)) % 2)
1953 {
1954 IdxRange[2][0] = nummodes[2] - 2;
1955 IdxRange[2][1] = 0;
1956 Incr[2] = -1;
1957 }
1958 else
1959 {
1960 IdxRange[2][0] = 1;
1961 IdxRange[2][1] = nummodes[2] - 1;
1962 Incr[2] = 1;
1963 }
1964 }
1965 else
1966 {
1967 IdxRange[2][0] = 2;
1968 IdxRange[2][1] = nummodes[2];
1969 Incr[2] = 1;
1970
1971 if (((int)(faceOrient - eDir1FwdDir1_Dir2FwdDir2)) % 2)
1972 {
1973 for (int i = 3; i < nummodes[2]; i += 2)
1974 {
1975 sign2[i] = -1;
1976 }
1977 }
1978 }
1979 }
1980 }
1981
1982 // Set the upper and lower bounds, and increment for the faces
1983 // involving the second coordinate direction.
1984 switch (fid)
1985 {
1986 case 1:
1987 {
1988 IdxRange[1][0] = 0;
1989 IdxRange[1][1] = 1;
1990 Incr[1] = 1;
1991 }
1992 break;
1993 case 3:
1994 {
1995 if (bType[1] == LibUtilities::eGLL_Lagrange)
1996 {
1997 IdxRange[1][0] = nummodes[1] - 1;
1998 IdxRange[1][1] = nummodes[1];
1999 Incr[1] = 1;
2000 }
2001 else
2002 {
2003 IdxRange[1][0] = 1;
2004 IdxRange[1][1] = 2;
2005 Incr[1] = 1;
2006 }
2007 }
2008 break;
2009 case 0:
2010 case 5:
2011 {
2012 if (bType[1] == LibUtilities::eGLL_Lagrange)
2013 {
2014 if (((int)(faceOrient - eDir1FwdDir1_Dir2FwdDir2)) % 2)
2015 {
2016 IdxRange[1][0] = nummodes[1] - 2;
2017 IdxRange[1][1] = 0;
2018 Incr[1] = -1;
2019 }
2020 else
2021 {
2022 IdxRange[1][0] = 1;
2023 IdxRange[1][1] = nummodes[1] - 1;
2024 Incr[1] = 1;
2025 }
2026 }
2027 else
2028 {
2029 IdxRange[1][0] = 2;
2030 IdxRange[1][1] = nummodes[1];
2031 Incr[1] = 1;
2032
2033 if (((int)(faceOrient - eDir1FwdDir1_Dir2FwdDir2)) % 2)
2034 {
2035 for (int i = 3; i < nummodes[1]; i += 2)
2036 {
2037 sign1[i] = -1;
2038 }
2039 }
2040 }
2041 }
2042 break;
2043 default: // case2: case4:
2044 {
2045 if (bType[1] == LibUtilities::eGLL_Lagrange)
2046 {
2047 if (((int)(faceOrient - eDir1FwdDir1_Dir2FwdDir2)) % 4 > 1)
2048 {
2049 IdxRange[1][0] = nummodes[1] - 2;
2050 IdxRange[1][1] = 0;
2051 Incr[1] = -1;
2052 }
2053 else
2054 {
2055 IdxRange[1][0] = 1;
2056 IdxRange[1][1] = nummodes[1] - 1;
2057 Incr[1] = 1;
2058 }
2059 }
2060 else
2061 {
2062 IdxRange[1][0] = 2;
2063 IdxRange[1][1] = nummodes[1];
2064 Incr[1] = 1;
2065
2066 if (((int)(faceOrient - eDir1FwdDir1_Dir2FwdDir2)) % 4 > 1)
2067 {
2068 for (int i = 3; i < nummodes[1]; i += 2)
2069 {
2070 sign1[i] = -1;
2071 }
2072 }
2073 }
2074 }
2075 }
2076
2077 switch (fid)
2078 {
2079 case 4:
2080 {
2081 IdxRange[0][0] = 0;
2082 IdxRange[0][1] = 1;
2083 Incr[0] = 1;
2084 }
2085 break;
2086 case 2:
2087 {
2088 if (bType[0] == LibUtilities::eGLL_Lagrange)
2089 {
2090 IdxRange[0][0] = nummodes[0] - 1;
2091 IdxRange[0][1] = nummodes[0];
2092 Incr[0] = 1;
2093 }
2094 else
2095 {
2096 IdxRange[0][0] = 1;
2097 IdxRange[0][1] = 2;
2098 Incr[0] = 1;
2099 }
2100 }
2101 break;
2102 default:
2103 {
2104 if (bType[0] == LibUtilities::eGLL_Lagrange)
2105 {
2106 if (((int)(faceOrient - eDir1FwdDir1_Dir2FwdDir2)) % 4 > 1)
2107 {
2108 IdxRange[0][0] = nummodes[0] - 2;
2109 IdxRange[0][1] = 0;
2110 Incr[0] = -1;
2111 }
2112 else
2113 {
2114 IdxRange[0][0] = 1;
2115 IdxRange[0][1] = nummodes[0] - 1;
2116 Incr[0] = 1;
2117 }
2118 }
2119 else
2120 {
2121 IdxRange[0][0] = 2;
2122 IdxRange[0][1] = nummodes[0];
2123 Incr[0] = 1;
2124
2125 if (((int)(faceOrient - eDir1FwdDir1_Dir2FwdDir2)) % 4 > 1)
2126 {
2127 for (int i = 3; i < nummodes[0]; i += 2)
2128 {
2129 sign0[i] = -1;
2130 }
2131 }
2132 }
2133 }
2134 }
2135
2136 int cnt = 0;
2137
2138 for (int r = IdxRange[2][0]; r != IdxRange[2][1]; r += Incr[2])
2139 {
2140 for (int q = IdxRange[1][0]; q != IdxRange[1][1]; q += Incr[1])
2141 {
2142 for (int p = IdxRange[0][0]; p != IdxRange[0][1]; p += Incr[0])
2143 {
2144 maparray[arrayindx[cnt]] =
2145 r * nummodes[0] * nummodes[1] + q * nummodes[0] + p;
2146 signarray[arrayindx[cnt++]] = sign0[p] * sign1[q] * sign2[r];
2147 }
2148 }
2149 }
2150}
int v_GetTraceIntNcoeffs(const int i) const override

References ASSERTL1, Nektar::StdRegions::eDir1FwdDir1_Dir2FwdDir2, Nektar::StdRegions::eDir1FwdDir2_Dir2FwdDir1, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::m_base, and v_GetTraceIntNcoeffs().

◆ v_GetTraceIntNcoeffs()

int Nektar::StdRegions::StdHexExp::v_GetTraceIntNcoeffs ( const int  i) const
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 660 of file StdHexExp.cpp.

661{
662 ASSERTL2((i >= 0) && (i <= 5), "face id is out of range");
663 if ((i == 0) || (i == 5))
664 {
665 return (GetBasisNumModes(0) - 2) * (GetBasisNumModes(1) - 2);
666 }
667 else if ((i == 1) || (i == 3))
668 {
669 return (GetBasisNumModes(0) - 2) * (GetBasisNumModes(2) - 2);
670 }
671 else
672 {
673 return (GetBasisNumModes(1) - 2) * (GetBasisNumModes(2) - 2);
674 }
675}

References ASSERTL2, and Nektar::StdRegions::StdExpansion::GetBasisNumModes().

Referenced by v_GetTraceInteriorToElementMap().

◆ v_GetTraceNcoeffs()

int Nektar::StdRegions::StdHexExp::v_GetTraceNcoeffs ( const int  i) const
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 643 of file StdHexExp.cpp.

644{
645 ASSERTL2((i >= 0) && (i <= 5), "face id is out of range");
646 if ((i == 0) || (i == 5))
647 {
648 return GetBasisNumModes(0) * GetBasisNumModes(1);
649 }
650 else if ((i == 1) || (i == 3))
651 {
652 return GetBasisNumModes(0) * GetBasisNumModes(2);
653 }
654 else
655 {
656 return GetBasisNumModes(1) * GetBasisNumModes(2);
657 }
658}

References ASSERTL2, and Nektar::StdRegions::StdExpansion::GetBasisNumModes().

◆ v_GetTraceNumModes()

void Nektar::StdRegions::StdHexExp::v_GetTraceNumModes ( const int  fid,
int &  numModes0,
int &  numModes1,
Orientation  faceOrient = eDir1FwdDir1_Dir2FwdDir2 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 779 of file StdHexExp.cpp.

781{
782 int nummodes[3] = {m_base[0]->GetNumModes(), m_base[1]->GetNumModes(),
783 m_base[2]->GetNumModes()};
784 switch (fid)
785 {
786 case 0:
787 case 5:
788 {
789 numModes0 = nummodes[0];
790 numModes1 = nummodes[1];
791 }
792 break;
793 case 1:
794 case 3:
795 {
796 numModes0 = nummodes[0];
797 numModes1 = nummodes[2];
798 }
799 break;
800 case 2:
801 case 4:
802 {
803 numModes0 = nummodes[1];
804 numModes1 = nummodes[2];
805 }
806 break;
807 default:
808 {
809 ASSERTL0(false, "fid out of range");
810 }
811 break;
812 }
813
814 if (faceOrient >= eDir1FwdDir2_Dir2FwdDir1)
815 {
816 std::swap(numModes0, numModes1);
817 }
818}

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

◆ v_GetTraceNumPoints()

int Nektar::StdRegions::StdHexExp::v_GetTraceNumPoints ( const int  i) const
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 677 of file StdHexExp.cpp.

678{
679 ASSERTL2(i >= 0 && i <= 5, "face id is out of range");
680
681 if (i == 0 || i == 5)
682 {
683 return m_base[0]->GetNumPoints() * m_base[1]->GetNumPoints();
684 }
685 else if (i == 1 || i == 3)
686 {
687 return m_base[0]->GetNumPoints() * m_base[2]->GetNumPoints();
688 }
689 else
690 {
691 return m_base[1]->GetNumPoints() * m_base[2]->GetNumPoints();
692 }
693}

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

◆ v_GetTracePointsKey()

LibUtilities::PointsKey Nektar::StdRegions::StdHexExp::v_GetTracePointsKey ( const int  i,
const int  j 
) const
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 695 of file StdHexExp.cpp.

697{
698 ASSERTL2(i >= 0 && i <= 5, "face id is out of range");
699 ASSERTL2(j == 0 || j == 1, "face direction is out of range");
700
701 if (i == 0 || i == 5)
702 {
703 return m_base[j]->GetPointsKey();
704 }
705 else if (i == 1 || i == 3)
706 {
707 return m_base[2 * j]->GetPointsKey();
708 }
709 else
710 {
711 return m_base[j + 1]->GetPointsKey();
712 }
713}

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

◆ v_GetVertexMap()

int Nektar::StdRegions::StdHexExp::v_GetVertexMap ( int  localVertexId,
bool  useCoeffPacking = false 
)
overrideprotectedvirtual

Expansions in each of the three dimensions must be of type LibUtilities::eModified_A or LibUtilities::eGLL_Lagrange.

Parameters
localVertexIdID of vertex (0..7)
Returns
Position of vertex in local numbering scheme.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 827 of file StdHexExp.cpp.

828{
831 "BasisType is not a boundary interior form");
834 "BasisType is not a boundary interior form");
837 "BasisType is not a boundary interior form");
838
839 ASSERTL1((localVertexId >= 0) && (localVertexId < 8),
840 "local vertex id must be between 0 and 7");
841
842 int p = 0;
843 int q = 0;
844 int r = 0;
845
846 // Retrieve the number of modes in each dimension.
847 int nummodes[3] = {m_base[0]->GetNumModes(), m_base[1]->GetNumModes(),
848 m_base[2]->GetNumModes()};
849
850 if (useCoeffPacking == true) // follow packing of coefficients i.e q,r,p
851 {
852 if (localVertexId > 3)
853 {
855 {
856 r = nummodes[2] - 1;
857 }
858 else
859 {
860 r = 1;
861 }
862 }
863
864 switch (localVertexId % 4)
865 {
866 case 0:
867 break;
868 case 1:
869 {
871 {
872 p = nummodes[0] - 1;
873 }
874 else
875 {
876 p = 1;
877 }
878 }
879 break;
880 case 2:
881 {
883 {
884 q = nummodes[1] - 1;
885 }
886 else
887 {
888 q = 1;
889 }
890 }
891 break;
892 case 3:
893 {
895 {
896 p = nummodes[0] - 1;
897 q = nummodes[1] - 1;
898 }
899 else
900 {
901 p = 1;
902 q = 1;
903 }
904 }
905 break;
906 }
907 }
908 else
909 {
910 // Right face (vertices 1,2,5,6)
911 if ((localVertexId % 4) % 3 > 0)
912 {
914 {
915 p = nummodes[0] - 1;
916 }
917 else
918 {
919 p = 1;
920 }
921 }
922 // Back face (vertices 2,3,6,7)
923 if (localVertexId % 4 > 1)
924 {
926 {
927 q = nummodes[1] - 1;
928 }
929 else
930 {
931 q = 1;
932 }
933 }
934
935 // Top face (vertices 4,5,6,7)
936 if (localVertexId > 3)
937 {
939 {
940 r = nummodes[2] - 1;
941 }
942 else
943 {
944 r = 1;
945 }
946 }
947 }
948 // Compute the local number.
949 return r * nummodes[0] * nummodes[1] + q * nummodes[0] + p;
950}

References ASSERTL1, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), and Nektar::StdRegions::StdExpansion::m_base.

◆ v_HelmholtzMatrixOp()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2293 of file StdHexExp.cpp.

2296{
2297 StdHexExp::v_HelmholtzMatrixOp_MatFree(inarray, outarray, mkey);
2298}
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)

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

◆ v_IProductWRTBase()

void Nektar::StdRegions::StdHexExp::v_IProductWRTBase ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

\( \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
inarray?
outarray?

Implements Nektar::StdRegions::StdExpansion.

Definition at line 314 of file StdHexExp.cpp.

316{
317 if (m_base[0]->Collocation() && m_base[1]->Collocation() &&
318 m_base[2]->Collocation())
319 {
320 MultiplyByQuadratureMetric(inarray, outarray);
321 }
322 else
323 {
324 StdHexExp::v_IProductWRTBase_SumFac(inarray, outarray);
325 }
326}
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multbyweights=true) override

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

◆ v_IProductWRTBase_SumFac()

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

Implementation of the sum-factorization inner product operation.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 331 of file StdHexExp.cpp.

334{
335 int nquad0 = m_base[0]->GetNumPoints();
336 int nquad1 = m_base[1]->GetNumPoints();
337 int nquad2 = m_base[2]->GetNumPoints();
338 int order0 = m_base[0]->GetNumModes();
339 int order1 = m_base[1]->GetNumModes();
340
341 Array<OneD, NekDouble> wsp(nquad0 * nquad1 * (nquad2 + order0) +
342 order0 * order1 * nquad2);
343
344 if (multiplybyweights)
345 {
346 Array<OneD, NekDouble> tmp(inarray.size());
347 MultiplyByQuadratureMetric(inarray, tmp);
348
350 m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
351 tmp, outarray, wsp, true, true, true);
352 }
353 else
354 {
356 m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
357 inarray, outarray, wsp, true, true, true);
358 }
359}
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_IProductWRTBase_SumFacKernel()

void Nektar::StdRegions::StdHexExp::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 
)
overrideprotectedvirtual

Implementation of the sum-factorisation inner product operation.

Todo:
Implement cases where only some directions are collocated.

Implements Nektar::StdRegions::StdExpansion3D.

Definition at line 365 of file StdHexExp.cpp.

372{
373 int nquad0 = m_base[0]->GetNumPoints();
374 int nquad1 = m_base[1]->GetNumPoints();
375 int nquad2 = m_base[2]->GetNumPoints();
376 int nmodes0 = m_base[0]->GetNumModes();
377 int nmodes1 = m_base[1]->GetNumModes();
378 int nmodes2 = m_base[2]->GetNumModes();
379
380 bool colldir0 = doCheckCollDir0 ? (m_base[0]->Collocation()) : false;
381 bool colldir1 = doCheckCollDir1 ? (m_base[1]->Collocation()) : false;
382 bool colldir2 = doCheckCollDir2 ? (m_base[2]->Collocation()) : false;
383
384 if (colldir0 && colldir1 && colldir2)
385 {
386 Vmath::Vcopy(m_ncoeffs, inarray.data(), 1, outarray.data(), 1);
387 }
388 else
389 {
390 ASSERTL1(wsp.size() >= nmodes0 * nquad2 * (nquad1 + nmodes1),
391 "Insufficient workspace size");
392
393 Array<OneD, NekDouble> tmp0 = wsp;
394 Array<OneD, NekDouble> tmp1 = wsp + nmodes0 * nquad1 * nquad2;
395
396 if (colldir0)
397 {
398 // reshuffle data for next operation.
399 for (int n = 0; n < nmodes0; ++n)
400 {
401 Vmath::Vcopy(nquad1 * nquad2, inarray.data() + n, nquad0,
402 tmp0.data() + nquad1 * nquad2 * n, 1);
403 }
404 }
405 else
406 {
407 Blas::Dgemm('T', 'N', nquad1 * nquad2, nmodes0, nquad0, 1.0,
408 inarray.data(), nquad0, base0.data(), nquad0, 0.0,
409 tmp0.data(), nquad1 * nquad2);
410 }
411
412 if (colldir1)
413 {
414 // reshuffle data for next operation.
415 for (int n = 0; n < nmodes1; ++n)
416 {
417 Vmath::Vcopy(nquad2 * nmodes0, tmp0.data() + n, nquad1,
418 tmp1.data() + nquad2 * nmodes0 * n, 1);
419 }
420 }
421 else
422 {
423 Blas::Dgemm('T', 'N', nquad2 * nmodes0, nmodes1, nquad1, 1.0,
424 tmp0.data(), nquad1, base1.data(), nquad1, 0.0,
425 tmp1.data(), nquad2 * nmodes0);
426 }
427
428 if (colldir2)
429 {
430 // reshuffle data for next operation.
431 for (int n = 0; n < nmodes2; ++n)
432 {
433 Vmath::Vcopy(nmodes0 * nmodes1, tmp1.data() + n, nquad2,
434 outarray.data() + nmodes0 * nmodes1 * n, 1);
435 }
436 }
437 else
438 {
439 Blas::Dgemm('T', 'N', nmodes0 * nmodes1, nmodes2, nquad2, 1.0,
440 tmp1.data(), nquad2, base2.data(), nquad2, 0.0,
441 outarray.data(), nmodes0 * nmodes1);
442 }
443 }
444}

References ASSERTL1, Blas::Dgemm(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, and Vmath::Vcopy().

◆ v_IProductWRTDerivBase()

void Nektar::StdRegions::StdHexExp::v_IProductWRTDerivBase ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 446 of file StdHexExp.cpp.

449{
450 StdHexExp::IProductWRTDerivBase_SumFac(dir, inarray, outarray);
451}
void IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)

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

◆ v_IProductWRTDerivBase_SumFac()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 453 of file StdHexExp.cpp.

456{
457 ASSERTL0((dir == 0) || (dir == 1) || (dir == 2),
458 "input dir is out of range");
459
460 int nquad1 = m_base[1]->GetNumPoints();
461 int nquad2 = m_base[2]->GetNumPoints();
462 int order0 = m_base[0]->GetNumModes();
463 int order1 = m_base[1]->GetNumModes();
464
465 // If outarray > inarray then no need for temporary storage.
466 Array<OneD, NekDouble> tmp = outarray;
467 if (outarray.size() < inarray.size())
468 {
469 tmp = Array<OneD, NekDouble>(inarray.size());
470 }
471
472 // Need workspace for sumfackernel though
473 Array<OneD, NekDouble> wsp(order0 * nquad2 * (nquad1 + order1));
474
475 // multiply by integration constants
476 MultiplyByQuadratureMetric(inarray, tmp);
477
478 // perform sum-factorisation
479 switch (dir)
480 {
481 case 0:
483 m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
484 m_base[2]->GetBdata(), tmp, outarray, wsp, false, true, true);
485 break;
486 case 1:
488 m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
489 m_base[2]->GetBdata(), tmp, outarray, wsp, true, false, true);
490 break;
491 case 2:
493 m_base[0]->GetBdata(), m_base[1]->GetBdata(),
494 m_base[2]->GetDbdata(), tmp, outarray, wsp, true, true, false);
495 break;
496 }
497}

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

◆ v_IsBoundaryInteriorExpansion()

bool Nektar::StdRegions::StdHexExp::v_IsBoundaryInteriorExpansion ( ) const
overrideprotectedvirtual

◆ v_LaplacianMatrixOp() [1/2]

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2270 of file StdHexExp.cpp.

2273{
2274 StdHexExp::v_LaplacianMatrixOp_MatFree(inarray, outarray, mkey);
2275}
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)

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

◆ v_LaplacianMatrixOp() [2/2]

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2277 of file StdHexExp.cpp.

2281{
2282 StdExpansion::LaplacianMatrixOp_MatFree(k1, k2, inarray, outarray, mkey);
2283}
void LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)

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

◆ v_LocCollapsedToLocCoord()

void Nektar::StdRegions::StdHexExp::v_LocCollapsedToLocCoord ( const Array< OneD, const NekDouble > &  eta,
Array< OneD, NekDouble > &  xi 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 507 of file StdHexExp.cpp.

509{
510 xi[0] = eta[0];
511 xi[1] = eta[1];
512 xi[2] = eta[2];
513}

◆ v_LocCoordToLocCollapsed()

void Nektar::StdRegions::StdHexExp::v_LocCoordToLocCollapsed ( const Array< OneD, const NekDouble > &  xi,
Array< OneD, NekDouble > &  eta 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 499 of file StdHexExp.cpp.

501{
502 eta[0] = xi[0];
503 eta[1] = xi[1];
504 eta[2] = xi[2];
505}

◆ v_MassMatrixOp()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2263 of file StdHexExp.cpp.

2266{
2267 StdExpansion::MassMatrixOp_MatFree(inarray, outarray, mkey);
2268}
void MassMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)

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

◆ v_MultiplyByStdQuadratureMetric()

void Nektar::StdRegions::StdHexExp::v_MultiplyByStdQuadratureMetric ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2300 of file StdHexExp.cpp.

2303{
2304 int nquad0 = m_base[0]->GetNumPoints();
2305 int nquad1 = m_base[1]->GetNumPoints();
2306 int nquad2 = m_base[2]->GetNumPoints();
2307 int nq01 = nquad0 * nquad1;
2308 int nq12 = nquad1 * nquad2;
2309
2310 const Array<OneD, const NekDouble> &w0 = m_base[0]->GetW();
2311 const Array<OneD, const NekDouble> &w1 = m_base[1]->GetW();
2312 const Array<OneD, const NekDouble> &w2 = m_base[2]->GetW();
2313
2314 for (int i = 0; i < nq12; ++i)
2315 {
2316 Vmath::Vmul(nquad0, inarray.data() + i * nquad0, 1, w0.data(), 1,
2317 outarray.data() + i * nquad0, 1);
2318 }
2319
2320 for (int i = 0; i < nq12; ++i)
2321 {
2322 Vmath::Smul(nquad0, w1[i % nquad1], outarray.data() + i * nquad0, 1,
2323 outarray.data() + i * nquad0, 1);
2324 }
2325
2326 for (int i = 0; i < nquad2; ++i)
2327 {
2328 Vmath::Smul(nq01, w2[i], outarray.data() + i * nq01, 1,
2329 outarray.data() + i * nq01, 1);
2330 }
2331}

References Nektar::StdRegions::StdExpansion::m_base, Vmath::Smul(), and Vmath::Vmul().

◆ v_NumBndryCoeffs()

int Nektar::StdRegions::StdHexExp::v_NumBndryCoeffs ( ) const
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 604 of file StdHexExp.cpp.

605{
608 "BasisType is not a boundary interior form");
611 "BasisType is not a boundary interior form");
614 "BasisType is not a boundary interior form");
615
616 int nmodes0 = m_base[0]->GetNumModes();
617 int nmodes1 = m_base[1]->GetNumModes();
618 int nmodes2 = m_base[2]->GetNumModes();
619
620 return (2 * (nmodes0 * nmodes1 + nmodes0 * nmodes2 + nmodes1 * nmodes2) -
621 4 * (nmodes0 + nmodes1 + nmodes2) + 8);
622}

References ASSERTL1, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), and Nektar::StdRegions::StdExpansion::m_base.

◆ v_NumDGBndryCoeffs()

int Nektar::StdRegions::StdHexExp::v_NumDGBndryCoeffs ( ) const
overrideprotectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 624 of file StdHexExp.cpp.

625{
628 "BasisType is not a boundary interior form");
631 "BasisType is not a boundary interior form");
634 "BasisType is not a boundary interior form");
635
636 int nmodes0 = m_base[0]->GetNumModes();
637 int nmodes1 = m_base[1]->GetNumModes();
638 int nmodes2 = m_base[2]->GetNumModes();
639
640 return 2 * (nmodes0 * nmodes1 + nmodes0 * nmodes2 + nmodes1 * nmodes2);
641}

References ASSERTL1, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), and Nektar::StdRegions::StdExpansion::m_base.

◆ v_PhysDeriv() [1/2]

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

Differentiation Methods.

For Hexahedral region can use the PhysTensorDeriv function defined under StdExpansion. Following tenserproduct:

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 73 of file StdHexExp.cpp.

77{
78 StdExpansion3D::PhysTensorDeriv(inarray, out_d0, out_d1, out_d2);
79}
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.

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

Referenced by v_StdPhysDeriv().

◆ v_PhysDeriv() [2/2]

void Nektar::StdRegions::StdHexExp::v_PhysDeriv ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual
Parameters
dirDirection in which to compute derivative. Valid values are 0, 1, 2.
inarrayInput array.
outarrayOutput array.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 87 of file StdHexExp.cpp.

90{
91 switch (dir)
92 {
93 case 0:
94 {
95 PhysDeriv(inarray, outarray, NullNekDouble1DArray,
97 }
98 break;
99 case 1:
100 {
101 PhysDeriv(inarray, NullNekDouble1DArray, outarray,
103 }
104 break;
105 case 2:
106 {
108 outarray);
109 }
110 break;
111 default:
112 {
113 ASSERTL1(false, "input dir is out of range");
114 }
115 break;
116 }
117}
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)
static Array< OneD, NekDouble > NullNekDouble1DArray

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

◆ v_PhysEvalFirstDeriv()

NekDouble Nektar::StdRegions::StdHexExp::v_PhysEvalFirstDeriv ( const Array< OneD, NekDouble > &  coord,
const Array< OneD, const NekDouble > &  inarray,
std::array< NekDouble, 3 > &  firstOrderDerivs 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1105 of file StdHexExp.cpp.

1109{
1110 return BaryTensorDeriv(coord, inarray, firstOrderDerivs);
1111}
NekDouble BaryTensorDeriv(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)

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

◆ v_PhysEvaluateBasis()

NekDouble Nektar::StdRegions::StdHexExp::v_PhysEvaluateBasis ( const Array< OneD, const NekDouble > &  coords,
int  mode 
)
finalprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 563 of file StdHexExp.cpp.

565{
566 ASSERTL2(coords[0] > -1 - NekConstants::kNekZeroTol, "coord[0] < -1");
567 ASSERTL2(coords[0] < 1 + NekConstants::kNekZeroTol, "coord[0] > 1");
568 ASSERTL2(coords[1] > -1 - NekConstants::kNekZeroTol, "coord[1] < -1");
569 ASSERTL2(coords[1] < 1 + NekConstants::kNekZeroTol, "coord[1] > 1");
570 ASSERTL2(coords[2] > -1 - NekConstants::kNekZeroTol, "coord[2] < -1");
571 ASSERTL2(coords[2] < 1 + NekConstants::kNekZeroTol, "coord[2] > 1");
572
573 const int nm0 = m_base[0]->GetNumModes();
574 const int nm1 = m_base[1]->GetNumModes();
575 const int mode2 = mode / (nm0 * nm1);
576 const int mode1 = (mode - mode2 * nm0 * nm1) / nm0;
577 const int mode0 = (mode - mode2 * nm0 * nm1) % nm0;
578
579 return StdExpansion::BaryEvaluateBasis<0>(coords[0], mode0) *
580 StdExpansion::BaryEvaluateBasis<1>(coords[1], mode1) *
581 StdExpansion::BaryEvaluateBasis<2>(coords[2], mode2);
582}
static const NekDouble kNekZeroTol

References ASSERTL2, Nektar::NekConstants::kNekZeroTol, and Nektar::StdRegions::StdExpansion::m_base.

◆ v_StdPhysDeriv()

void Nektar::StdRegions::StdHexExp::v_StdPhysDeriv ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d0,
Array< OneD, NekDouble > &  out_d1,
Array< OneD, NekDouble > &  out_d2 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 119 of file StdHexExp.cpp.

123{
124 StdHexExp::v_PhysDeriv(inarray, out_d0, out_d1, out_d2);
125}
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.
Definition StdHexExp.cpp:73

References v_PhysDeriv().

◆ v_SVVLaplacianFilter()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2333 of file StdHexExp.cpp.

2335{
2336 // Generate an orthonogal expansion
2337 int qa = m_base[0]->GetNumPoints();
2338 int qb = m_base[1]->GetNumPoints();
2339 int qc = m_base[2]->GetNumPoints();
2340 int nmodes_a = m_base[0]->GetNumModes();
2341 int nmodes_b = m_base[1]->GetNumModes();
2342 int nmodes_c = m_base[2]->GetNumModes();
2343 // Declare orthogonal basis.
2344 LibUtilities::PointsKey pa(qa, m_base[0]->GetPointsType());
2345 LibUtilities::PointsKey pb(qb, m_base[1]->GetPointsType());
2346 LibUtilities::PointsKey pc(qc, m_base[2]->GetPointsType());
2347
2348 LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A, nmodes_a, pa);
2349 LibUtilities::BasisKey Bb(LibUtilities::eOrtho_A, nmodes_b, pb);
2350 LibUtilities::BasisKey Bc(LibUtilities::eOrtho_A, nmodes_c, pc);
2351 StdHexExp OrthoExp(Ba, Bb, Bc);
2352
2353 Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
2354 int cnt = 0;
2355
2356 // project onto modal space.
2357 OrthoExp.FwdTrans(array, orthocoeffs);
2358
2359 if (mkey.ConstFactorExists(eFactorSVVPowerKerDiffCoeff))
2360 {
2361 // Rodrigo's power kernel
2362 NekDouble cutoff = mkey.GetConstFactor(eFactorSVVCutoffRatio);
2363 NekDouble SvvDiffCoeff =
2364 mkey.GetConstFactor(eFactorSVVPowerKerDiffCoeff) *
2365 mkey.GetConstFactor(eFactorSVVDiffCoeff);
2366
2367 for (int i = 0; i < nmodes_a; ++i)
2368 {
2369 for (int j = 0; j < nmodes_b; ++j)
2370 {
2371 NekDouble fac1 = std::max(
2372 pow((1.0 * i) / (nmodes_a - 1), cutoff * nmodes_a),
2373 pow((1.0 * j) / (nmodes_b - 1), cutoff * nmodes_b));
2374
2375 for (int k = 0; k < nmodes_c; ++k)
2376 {
2377 NekDouble fac =
2378 std::max(fac1, pow((1.0 * k) / (nmodes_c - 1),
2379 cutoff * nmodes_c));
2380
2381 orthocoeffs[cnt] *= SvvDiffCoeff * fac;
2382 cnt++;
2383 }
2384 }
2385 }
2386 }
2387 else if (mkey.ConstFactorExists(
2388 eFactorSVVDGKerDiffCoeff)) // Rodrigo/Mansoor's DG Kernel
2389 {
2390 NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDGKerDiffCoeff) *
2391 mkey.GetConstFactor(eFactorSVVDiffCoeff);
2392
2393 int max_abc = max(nmodes_a - kSVVDGFiltermodesmin,
2394 nmodes_b - kSVVDGFiltermodesmin);
2395 max_abc = max(max_abc, nmodes_c - kSVVDGFiltermodesmin);
2396 // clamp max_abc
2397 max_abc = max(max_abc, 0);
2398 max_abc = min(max_abc, kSVVDGFiltermodesmax - kSVVDGFiltermodesmin);
2399
2400 for (int i = 0; i < nmodes_a; ++i)
2401 {
2402 for (int j = 0; j < nmodes_b; ++j)
2403 {
2404 int maxij = max(i, j);
2405
2406 for (int k = 0; k < nmodes_c; ++k)
2407 {
2408 int maxijk = max(maxij, k);
2409 maxijk = min(maxijk, kSVVDGFiltermodesmax - 1);
2410
2411 orthocoeffs[cnt] *=
2412 SvvDiffCoeff * kSVVDGFilter[max_abc][maxijk];
2413 cnt++;
2414 }
2415 }
2416 }
2417 }
2418 else
2419 {
2420
2421 int cutoff = (int)(mkey.GetConstFactor(eFactorSVVCutoffRatio) *
2422 min(nmodes_a, nmodes_b));
2423 NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDiffCoeff);
2424 // Filter just trilinear space
2425 int nmodes = max(nmodes_a, nmodes_b);
2426 nmodes = max(nmodes, nmodes_c);
2427
2428 Array<OneD, NekDouble> fac(nmodes, 1.0);
2429 for (int j = cutoff; j < nmodes; ++j)
2430 {
2431 fac[j] = fabs((j - nmodes) / ((NekDouble)(j - cutoff + 1.0)));
2432 fac[j] *= fac[j]; // added this line to conform with equation
2433 }
2434
2435 for (int i = 0; i < nmodes_a; ++i)
2436 {
2437 for (int j = 0; j < nmodes_b; ++j)
2438 {
2439 for (int k = 0; k < nmodes_c; ++k)
2440 {
2441 if ((i >= cutoff) || (j >= cutoff) || (k >= cutoff))
2442 {
2443 orthocoeffs[i * nmodes_a * nmodes_b + j * nmodes_c +
2444 k] *=
2445 (SvvDiffCoeff * exp(-(fac[i] + fac[j] + fac[k])));
2446 }
2447 else
2448 {
2449 orthocoeffs[i * nmodes_a * nmodes_b + j * nmodes_c +
2450 k] *= 0.0;
2451 }
2452 }
2453 }
2454 }
2455 }
2456
2457 // backward transform to physical space
2458 OrthoExp.BwdTrans(orthocoeffs, array);
2459}
const int kSVVDGFiltermodesmin
const int kSVVDGFiltermodesmax
const NekDouble kSVVDGFilter[9][11]
scalarT< T > min(scalarT< T > lhs, scalarT< T > rhs)
Definition scalar.hpp:300

References Nektar::StdRegions::StdExpansion::BwdTrans(), Nektar::StdRegions::StdMatrixKey::ConstFactorExists(), Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::eFactorSVVDGKerDiffCoeff, Nektar::StdRegions::eFactorSVVDiffCoeff, Nektar::StdRegions::eFactorSVVPowerKerDiffCoeff, Nektar::LibUtilities::eOrtho_A, Nektar::StdRegions::StdExpansion::FwdTrans(), Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::StdRegions::kSVVDGFilter, Nektar::StdRegions::kSVVDGFiltermodesmax, Nektar::StdRegions::kSVVDGFiltermodesmin, Nektar::StdRegions::StdExpansion::m_base, tinysimd::max(), and tinysimd::min().

◆ v_WeakDerivMatrixOp()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2285 of file StdHexExp.cpp.

2289{
2290 StdExpansion::WeakDerivMatrixOp_MatFree(i, inarray, outarray, mkey);
2291}
void WeakDerivMatrixOp_MatFree(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)

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