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

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. More...
 
virtual void v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
virtual void v_StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2) override
 
virtual void v_StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2) override
 
void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multbyweights=true) override
 
void v_IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2) override
 
void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
virtual void v_LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta) override
 
virtual void v_LocCollapsedToLocCoord (const Array< OneD, const NekDouble > &eta, Array< OneD, NekDouble > &xi) override
 
virtual void v_FillMode (const int mode, Array< OneD, NekDouble > &outarray) override
 
NekDouble v_PhysEvaluateBasis (const Array< OneD, const NekDouble > &coords, int mode) final override
 
virtual NekDouble v_PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
 
virtual int v_GetNverts () const override
 
virtual int v_GetNedges () const override
 
virtual int v_GetNtraces () const override
 
virtual LibUtilities::ShapeType v_DetShapeType () const override
 
virtual int v_NumBndryCoeffs () const override
 
virtual int v_NumDGBndryCoeffs () const override
 
virtual int v_GetTraceNcoeffs (const int i) const override
 
virtual int v_GetTraceIntNcoeffs (const int i) const override
 
virtual int v_GetTraceNumPoints (const int i) const override
 
virtual LibUtilities::PointsKey v_GetTracePointsKey (const int i, const int j) const override
 
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset) override
 
virtual const LibUtilities::BasisKey v_GetTraceBasisKey (const int i, const int k) const override
 
virtual bool v_IsBoundaryInteriorExpansion () const override
 
virtual void v_GetCoords (Array< OneD, NekDouble > &coords_x, Array< OneD, NekDouble > &coords_y, Array< OneD, NekDouble > &coords_z) override
 
virtual void v_GetTraceNumModes (const int fid, int &numModes0, int &numModes1, Orientation faceOrient=eDir1FwdDir1_Dir2FwdDir2) override
 
virtual int v_GetEdgeNcoeffs (const int i) const override
 
virtual int v_GetVertexMap (int localVertexId, bool useCoeffPacking=false) override
 
virtual void v_GetInteriorMap (Array< OneD, unsigned int > &outarray) override
 
virtual void v_GetBoundaryMap (Array< OneD, unsigned int > &outarray) override
 
virtual void v_GetTraceCoeffMap (const unsigned int fid, Array< OneD, unsigned int > &maparray) override
 
virtual void v_GetElmtTraceToTraceMap (const unsigned int fid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation faceOrient, int P, int Q) override
 
virtual 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
 
virtual DNekMatSharedPtr v_GenMatrix (const StdMatrixKey &mkey) override
 
virtual DNekMatSharedPtr v_CreateStdMatrix (const StdMatrixKey &mkey) override
 
virtual 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
 
virtual 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
virtual NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals) override
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
virtual NekDouble v_PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals) override
 
virtual NekDouble v_PhysEvaluate (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs) override
 
virtual void v_BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)=0
 
virtual void v_IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)=0
 
virtual void v_LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
virtual void v_HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override
 
virtual NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray) override
 Integrates the specified function over the domain. More...
 
virtual int v_GetNedges (void) const
 
virtual int v_GetEdgeNcoeffs (const int i) const
 
NekDouble BaryTensorDeriv (const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
 
virtual void v_GetEdgeInteriorToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient=eForwards)
 
virtual void v_GetTraceToElementMap (const int tid, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, Orientation traceOrient, int P, int Q) override
 
virtual void v_GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat) override
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion
DNekMatSharedPtr CreateStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr CreateStdStaticCondMatrix (const StdMatrixKey &mkey)
 Create the static condensation of a matrix when using a boundary interior decomposition. More...
 
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDirectionalDerivBase_SumFac (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GeneralMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
 
void LaplacianMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp_MatFree (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void HelmholtzMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
 
virtual NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
virtual void v_GenStdMatBwdDeriv (const int dir, DNekMatSharedPtr &mat)
 
virtual void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals, NekDouble &deriv, NekDouble &deriv2)
 This function performs the barycentric interpolation of the polynomial stored in coord at a point physvals using barycentric interpolation weights in direction. More...
 
template<int DIR>
NekDouble BaryEvaluateBasis (const NekDouble &coord, const int &mode)
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals)
 Helper function to pass an unused value by reference into BaryEvaluate. More...
 
template<int DIR, bool DERIV = false, bool DERIV2 = false>
NekDouble BaryEvaluate (const NekDouble &coord, const NekDouble *physvals, NekDouble &deriv)
 

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 47 of file StdHexExp.h.

Constructor & Destructor Documentation

◆ StdHexExp() [1/3]

Nektar::StdRegions::StdHexExp::StdHexExp ( )

Definition at line 47 of file StdHexExp.cpp.

48{
49}

◆ StdHexExp() [2/3]

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

Definition at line 51 of file StdHexExp.cpp.

54 : StdExpansion(Ba.GetNumModes() * Bb.GetNumModes() * Bc.GetNumModes(), 3,
55 Ba, Bb, Bc),
56 StdExpansion3D(Ba.GetNumModes() * Bb.GetNumModes() * Bc.GetNumModes(), Ba,
57 Bb, Bc)
58{
59}
StdExpansion()
Default Constructor.

◆ StdHexExp() [3/3]

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

Definition at line 61 of file StdHexExp.cpp.

62{
63}

◆ ~StdHexExp()

Nektar::StdRegions::StdHexExp::~StdHexExp ( )
overridevirtual

Definition at line 65 of file StdHexExp.cpp.

66{
67}

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 167 of file StdHexExp.cpp.

169{
172 "Basis[1] is not a general tensor type");
173
176 "Basis[2] is not a general tensor type");
177
178 if (m_base[0]->Collocation() && m_base[1]->Collocation() &&
179 m_base[2]->Collocation())
180 {
182 m_base[2]->GetNumPoints(),
183 inarray, 1, outarray, 1);
184 }
185 else
186 {
187 StdHexExp::BwdTrans_SumFac(inarray, outarray);
188 }
189}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:162
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.
Definition: StdExpansion.h:224
Array< OneD, LibUtilities::BasisSharedPtr > m_base
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:51
@ eModified_C
Principle Modified Functions .
Definition: BasisType.h:52
@ eOrtho_C
Principle Orthogonal Functions .
Definition: BasisType.h:48
@ eOrtho_B
Principle Orthogonal Functions .
Definition: BasisType.h:46
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191

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 194 of file StdHexExp.cpp.

196{
197 Array<OneD, NekDouble> wsp(
198 m_base[0]->GetNumPoints() * m_base[2]->GetNumModes() *
199 (m_base[1]->GetNumModes() + m_base[1]->GetNumPoints())); // FIX THIS
200
201 BwdTrans_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetBdata(),
202 m_base[2]->GetBdata(), inarray, outarray, wsp, true,
203 true, true);
204}
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 219 of file StdHexExp.cpp.

226{
227 int nquad0 = m_base[0]->GetNumPoints();
228 int nquad1 = m_base[1]->GetNumPoints();
229 int nquad2 = m_base[2]->GetNumPoints();
230 int nmodes0 = m_base[0]->GetNumModes();
231 int nmodes1 = m_base[1]->GetNumModes();
232 int nmodes2 = m_base[2]->GetNumModes();
233
234 // Check if using collocation, if requested.
235 bool colldir0 = doCheckCollDir0 ? (m_base[0]->Collocation()) : false;
236 bool colldir1 = doCheckCollDir1 ? (m_base[1]->Collocation()) : false;
237 bool colldir2 = doCheckCollDir2 ? (m_base[2]->Collocation()) : false;
238
239 // If collocation in all directions, Physical values at quadrature
240 // points is just a copy of the modes.
241 if (colldir0 && colldir1 && colldir2)
242 {
243 Vmath::Vcopy(m_ncoeffs, inarray.get(), 1, outarray.get(), 1);
244 }
245 else
246 {
247 // Check sufficiently large workspace.
248 ASSERTL1(wsp.size() >= nquad0 * nmodes2 * (nmodes1 + nquad1),
249 "Workspace size is not sufficient");
250
251 // Assign second half of workspace for 2nd DGEMM operation.
252 Array<OneD, NekDouble> wsp2 = wsp + nquad0 * nmodes1 * nmodes2;
253
254 // BwdTrans in each direction using DGEMM
255 Blas::Dgemm('T', 'T', nmodes1 * nmodes2, nquad0, nmodes0, 1.0,
256 &inarray[0], nmodes0, base0.get(), nquad0, 0.0, &wsp[0],
257 nmodes1 * nmodes2);
258 Blas::Dgemm('T', 'T', nquad0 * nmodes2, nquad1, nmodes1, 1.0, &wsp[0],
259 nmodes1, base1.get(), nquad1, 0.0, &wsp2[0],
260 nquad0 * nmodes2);
261 Blas::Dgemm('T', 'T', nquad0 * nquad1, nquad2, nmodes2, 1.0, &wsp2[0],
262 nmodes2, base2.get(), nquad2, 0.0, &outarray[0],
263 nquad0 * nquad1);
264 }
265}
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:385

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 735 of file StdHexExp.cpp.

737{
738 int nmodes = nummodes[modes_offset] * nummodes[modes_offset + 1] *
739 nummodes[modes_offset + 2];
740 modes_offset += 3;
741
742 return nmodes;
743}

◆ v_CreateStdMatrix()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 2279 of file StdHexExp.cpp.

2280{
2281 return v_GenMatrix(mkey);
2282}
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override
Definition: StdHexExp.cpp:2191

References v_GenMatrix().

◆ v_DetShapeType()

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

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 619 of file StdHexExp.cpp.

620{
622}

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 2482 of file StdHexExp.cpp.

2486{
2487 // Generate an orthogonal expansion
2488 int qa = m_base[0]->GetNumPoints();
2489 int qb = m_base[1]->GetNumPoints();
2490 int qc = m_base[2]->GetNumPoints();
2491 int nmodesA = m_base[0]->GetNumModes();
2492 int nmodesB = m_base[1]->GetNumModes();
2493 int nmodesC = m_base[2]->GetNumModes();
2494 int P = nmodesA - 1;
2495 int Q = nmodesB - 1;
2496 int R = nmodesC - 1;
2497
2498 // Declare orthogonal basis.
2499 LibUtilities::PointsKey pa(qa, m_base[0]->GetPointsType());
2500 LibUtilities::PointsKey pb(qb, m_base[1]->GetPointsType());
2501 LibUtilities::PointsKey pc(qc, m_base[2]->GetPointsType());
2502
2503 LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A, nmodesA, pa);
2504 LibUtilities::BasisKey Bb(LibUtilities::eOrtho_A, nmodesB, pb);
2505 LibUtilities::BasisKey Bc(LibUtilities::eOrtho_A, nmodesC, pc);
2506 StdHexExp OrthoExp(Ba, Bb, Bc);
2507
2508 // Cutoff
2509 int Pcut = cutoff * P;
2510 int Qcut = cutoff * Q;
2511 int Rcut = cutoff * R;
2512
2513 // Project onto orthogonal space.
2514 Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
2515 OrthoExp.FwdTrans(array, orthocoeffs);
2516
2517 //
2518 NekDouble fac, fac1, fac2, fac3;
2519 int index = 0;
2520 for (int i = 0; i < nmodesA; ++i)
2521 {
2522 for (int j = 0; j < nmodesB; ++j)
2523 {
2524 for (int k = 0; k < nmodesC; ++k, ++index)
2525 {
2526 // to filter out only the "high-modes"
2527 if (i > Pcut || j > Qcut || k > Rcut)
2528 {
2529 fac1 = (NekDouble)(i - Pcut) / ((NekDouble)(P - Pcut));
2530 fac2 = (NekDouble)(j - Qcut) / ((NekDouble)(Q - Qcut));
2531 fac3 = (NekDouble)(k - Rcut) / ((NekDouble)(R - Rcut));
2532 fac = max(max(fac1, fac2), fac3);
2533 fac = pow(fac, exponent);
2534 orthocoeffs[index] *= exp(-alpha * fac);
2535 }
2536 }
2537 }
2538 }
2539
2540 // backward transform to physical space
2541 OrthoExp.BwdTrans(orthocoeffs, array);
2542}
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:211
@ P
Monomial polynomials .
Definition: BasisType.h:64
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:44
double NekDouble

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, 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 538 of file StdHexExp.cpp.

539{
540 int nquad0 = m_base[0]->GetNumPoints();
541 int nquad1 = m_base[1]->GetNumPoints();
542 int nquad2 = m_base[2]->GetNumPoints();
543
544 Array<OneD, const NekDouble> base0 = m_base[0]->GetBdata();
545 Array<OneD, const NekDouble> base1 = m_base[1]->GetBdata();
546 Array<OneD, const NekDouble> base2 = m_base[2]->GetBdata();
547
548 int btmp0 = m_base[0]->GetNumModes();
549 int btmp1 = m_base[1]->GetNumModes();
550 int mode2 = mode / (btmp0 * btmp1);
551 int mode1 = (mode - mode2 * btmp0 * btmp1) / btmp0;
552 int mode0 = (mode - mode2 * btmp0 * btmp1) % btmp0;
553
554 ASSERTL2(mode == mode2 * btmp0 * btmp1 + mode1 * btmp0 + mode0,
555 "Mode lookup failed.");
556 ASSERTL2(mode < m_ncoeffs,
557 "Calling argument mode is larger than total expansion "
558 "order");
559
560 for (int i = 0; i < nquad1 * nquad2; ++i)
561 {
562 Vmath::Vcopy(nquad0, (NekDouble *)(base0.get() + mode0 * nquad0), 1,
563 &outarray[0] + i * nquad0, 1);
564 }
565
566 for (int j = 0; j < nquad2; ++j)
567 {
568 for (int i = 0; i < nquad0; ++i)
569 {
570 Vmath::Vmul(nquad1, (NekDouble *)(base1.get() + mode1 * nquad1), 1,
571 &outarray[0] + i + j * nquad0 * nquad1, nquad0,
572 &outarray[0] + i + j * nquad0 * nquad1, nquad0);
573 }
574 }
575
576 for (int i = 0; i < nquad2; i++)
577 {
578 Blas::Dscal(nquad0 * nquad1, base2[mode2 * nquad2 + i],
579 &outarray[0] + i * nquad0 * nquad1, 1);
580 }
581}
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:272
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition: Blas.hpp:151
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:207

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.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 276 of file StdHexExp.cpp.

278{
279 // If using collocation expansion, coefficients match physical
280 // data points so just do a direct copy.
281 if ((m_base[0]->Collocation()) && (m_base[1]->Collocation()) &&
282 (m_base[2]->Collocation()))
283 {
284 Vmath::Vcopy(GetNcoeffs(), &inarray[0], 1, &outarray[0], 1);
285 }
286 else
287 {
288 // Compute B^TWu
289 IProductWRTBase(inarray, outarray);
290
291 // get Mass matrix inverse
292 StdMatrixKey masskey(eInvMass, DetShapeType(), *this);
293 DNekMatSharedPtr matsys = GetStdMatrix(masskey);
294
295 // copy inarray in case inarray == outarray
296 DNekVec in(m_ncoeffs, outarray);
297 DNekVec out(m_ncoeffs, outarray, eWrapper);
298
299 // Solve for coefficients.
300 out = (*matsys) * in;
301 }
302}
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:130
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:609
void IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
this function calculates the inner product of a given function f with the different modes of the expa...
Definition: StdExpansion.h:534
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:373
NekVector< NekDouble > DNekVec
Definition: NekTypeDefs.hpp:48
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75

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.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 2191 of file StdHexExp.cpp.

2192{
2193
2194 MatrixType mtype = mkey.GetMatrixType();
2195
2196 DNekMatSharedPtr Mat;
2197
2198 switch (mtype)
2199 {
2201 {
2202 int nq0 = m_base[0]->GetNumPoints();
2203 int nq1 = m_base[1]->GetNumPoints();
2204 int nq2 = m_base[2]->GetNumPoints();
2205 int nq;
2206
2207 // take definition from key
2208 if (mkey.ConstFactorExists(eFactorConst))
2209 {
2210 nq = (int)mkey.GetConstFactor(eFactorConst);
2211 }
2212 else
2213 {
2214 nq = max(nq0, max(nq1, nq2));
2215 }
2216
2217 int neq =
2219 Array<OneD, Array<OneD, NekDouble>> coords(neq);
2220 Array<OneD, NekDouble> coll(3);
2221 Array<OneD, DNekMatSharedPtr> I(3);
2222 Array<OneD, NekDouble> tmp(nq0);
2223
2224 Mat =
2225 MemoryManager<DNekMat>::AllocateSharedPtr(neq, nq0 * nq1 * nq2);
2226 int cnt = 0;
2227
2228 for (int i = 0; i < nq; ++i)
2229 {
2230 for (int j = 0; j < nq; ++j)
2231 {
2232 for (int k = 0; k < nq; ++k, ++cnt)
2233 {
2234 coords[cnt] = Array<OneD, NekDouble>(3);
2235 coords[cnt][0] = -1.0 + 2 * k / (NekDouble)(nq - 1);
2236 coords[cnt][1] = -1.0 + 2 * j / (NekDouble)(nq - 1);
2237 coords[cnt][2] = -1.0 + 2 * i / (NekDouble)(nq - 1);
2238 }
2239 }
2240 }
2241
2242 for (int i = 0; i < neq; ++i)
2243 {
2244 LocCoordToLocCollapsed(coords[i], coll);
2245
2246 I[0] = m_base[0]->GetI(coll);
2247 I[1] = m_base[1]->GetI(coll + 1);
2248 I[2] = m_base[2]->GetI(coll + 2);
2249
2250 // interpolate first coordinate direction
2251 NekDouble fac;
2252 for (int k = 0; k < nq2; ++k)
2253 {
2254 for (int j = 0; j < nq1; ++j)
2255 {
2256
2257 fac = (I[1]->GetPtr())[j] * (I[2]->GetPtr())[k];
2258 Vmath::Smul(nq0, fac, I[0]->GetPtr(), 1, tmp, 1);
2259
2260 Vmath::Vcopy(nq0, &tmp[0], 1,
2261 Mat->GetRawPtr() + k * nq0 * nq1 * neq +
2262 j * nq0 * neq + i,
2263 neq);
2264 }
2265 }
2266 }
2267 }
2268 break;
2269 default:
2270 {
2272 }
2273 break;
2274 }
2275
2276 return Mat;
2277}
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)
Definition: ShapeType.hpp:158
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:245

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, 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 1038 of file StdHexExp.cpp.

1039{
1042 "BasisType is not a boundary interior form");
1045 "BasisType is not a boundary interior form");
1048 "BasisType is not a boundary interior form");
1049
1050 int i;
1051 int nummodes[3] = {m_base[0]->GetNumModes(), m_base[1]->GetNumModes(),
1052 m_base[2]->GetNumModes()};
1053
1054 int nBndCoeffs = NumBndryCoeffs();
1055
1056 if (outarray.size() != nBndCoeffs)
1057 {
1058 outarray = Array<OneD, unsigned int>(nBndCoeffs);
1059 }
1060
1061 const LibUtilities::BasisType Btype[3] = {GetBasisType(0), GetBasisType(1),
1062 GetBasisType(2)};
1063
1064 int p, q, r;
1065 int cnt = 0;
1066
1067 int BndIdx[3][2];
1068 int IntIdx[3][2];
1069
1070 for (i = 0; i < 3; i++)
1071 {
1072 BndIdx[i][0] = 0;
1073
1074 if (Btype[i] == LibUtilities::eModified_A)
1075 {
1076 BndIdx[i][1] = 1;
1077 IntIdx[i][0] = 2;
1078 IntIdx[i][1] = nummodes[i];
1079 }
1080 else
1081 {
1082 BndIdx[i][1] = nummodes[i] - 1;
1083 IntIdx[i][0] = 1;
1084 IntIdx[i][1] = nummodes[i] - 1;
1085 }
1086 }
1087
1088 for (i = 0; i < 2; i++)
1089 {
1090 r = BndIdx[2][i];
1091 for (q = 0; q < nummodes[1]; q++)
1092 {
1093 for (p = 0; p < nummodes[0]; p++)
1094 {
1095 outarray[cnt++] =
1096 r * nummodes[0] * nummodes[1] + q * nummodes[0] + p;
1097 }
1098 }
1099 }
1100
1101 for (r = IntIdx[2][0]; r < IntIdx[2][1]; r++)
1102 {
1103 for (i = 0; i < 2; i++)
1104 {
1105 q = BndIdx[1][i];
1106 for (p = 0; p < nummodes[0]; p++)
1107 {
1108 outarray[cnt++] =
1109 r * nummodes[0] * nummodes[1] + q * nummodes[0] + p;
1110 }
1111 }
1112
1113 for (q = IntIdx[1][0]; q < IntIdx[1][1]; q++)
1114 {
1115 for (i = 0; i < 2; i++)
1116 {
1117 p = BndIdx[0][i];
1118 outarray[cnt++] =
1119 r * nummodes[0] * nummodes[1] + q * nummodes[0] + p;
1120 }
1121 }
1122 }
1123
1124 sort(outarray.get(), outarray.get() + nBndCoeffs);
1125}
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:58
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:50
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, Nektar::StdRegions::StdExpansion::NumBndryCoeffs(), CellMLToNektar.cellml_metadata::p, and Nektar::UnitTests::q().

◆ 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.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 773 of file StdHexExp.cpp.

776{
777 Array<OneD, const NekDouble> eta_x = m_base[0]->GetZ();
778 Array<OneD, const NekDouble> eta_y = m_base[1]->GetZ();
779 Array<OneD, const NekDouble> eta_z = m_base[2]->GetZ();
780 int Qx = GetNumPoints(0);
781 int Qy = GetNumPoints(1);
782 int Qz = GetNumPoints(2);
783
784 // Convert collapsed coordinates into cartesian coordinates:
785 // eta --> xi
786 for (int k = 0; k < Qz; ++k)
787 {
788 for (int j = 0; j < Qy; ++j)
789 {
790 for (int i = 0; i < Qx; ++i)
791 {
792 int s = i + Qx * (j + Qy * k);
793 xi_x[s] = eta_x[i];
794 xi_y[s] = eta_y[j];
795 xi_z[s] = eta_z[k];
796 }
797 }
798 }
799}

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 1556 of file StdHexExp.cpp.

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

◆ v_GetEdgeNcoeffs()

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

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 2173 of file StdHexExp.cpp.

2174{
2175 ASSERTL2((i >= 0) && (i <= 11), "edge id is out of range");
2176
2177 if ((i == 0) || (i == 2) || (i == 8) || (i == 10))
2178 {
2179 return GetBasisNumModes(0);
2180 }
2181 else if ((i == 1) || (i == 3) || (i == 9) || (i == 11))
2182 {
2183 return GetBasisNumModes(1);
2184 }
2185 else
2186 {
2187 return GetBasisNumModes(2);
2188 }
2189}
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:175

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 1261 of file StdHexExp.cpp.

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

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 977 of file StdHexExp.cpp.

978{
981 "BasisType is not a boundary interior form");
984 "BasisType is not a boundary interior form");
987 "BasisType is not a boundary interior form");
988
989 int i;
990 int nummodes[3] = {m_base[0]->GetNumModes(), m_base[1]->GetNumModes(),
991 m_base[2]->GetNumModes()};
992
993 int nIntCoeffs = m_ncoeffs - NumBndryCoeffs();
994
995 if (outarray.size() != nIntCoeffs)
996 {
997 outarray = Array<OneD, unsigned int>(nIntCoeffs);
998 }
999
1000 const LibUtilities::BasisType Btype[3] = {GetBasisType(0), GetBasisType(1),
1001 GetBasisType(2)};
1002
1003 int p, q, r;
1004 int cnt = 0;
1005
1006 int IntIdx[3][2];
1007
1008 for (i = 0; i < 3; i++)
1009 {
1010 if (Btype[i] == LibUtilities::eModified_A)
1011 {
1012 IntIdx[i][0] = 2;
1013 IntIdx[i][1] = nummodes[i];
1014 }
1015 else
1016 {
1017 IntIdx[i][0] = 1;
1018 IntIdx[i][1] = nummodes[i] - 1;
1019 }
1020 }
1021
1022 for (r = IntIdx[2][0]; r < IntIdx[2][1]; r++)
1023 {
1024 for (q = IntIdx[1][0]; q < IntIdx[1][1]; q++)
1025 {
1026 for (p = IntIdx[0][0]; p < IntIdx[0][1]; p++)
1027 {
1028 outarray[cnt++] =
1029 r * nummodes[0] * nummodes[1] + q * nummodes[0] + p;
1030 }
1031 }
1032 }
1033}

References ASSERTL1, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::StdRegions::StdExpansion::NumBndryCoeffs(), CellMLToNektar.cellml_metadata::p, and Nektar::UnitTests::q().

◆ v_GetNedges()

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

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 609 of file StdHexExp.cpp.

610{
611 return 12;
612}

◆ v_GetNtraces()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 614 of file StdHexExp.cpp.

615{
616 return 6;
617}

◆ v_GetNverts()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 604 of file StdHexExp.cpp.

605{
606 return 8;
607}

◆ v_GetSimplexEquiSpacedConnectivity()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2544 of file StdHexExp.cpp.

2546{
2547 boost::ignore_unused(standard);
2548
2549 int np0 = m_base[0]->GetNumPoints();
2550 int np1 = m_base[1]->GetNumPoints();
2551 int np2 = m_base[2]->GetNumPoints();
2552 int np = max(np0, max(np1, np2));
2553
2554 conn = Array<OneD, int>(6 * (np - 1) * (np - 1) * (np - 1));
2555
2556 int row = 0;
2557 int rowp1 = 0;
2558 int cnt = 0;
2559 int plane = 0;
2560 for (int i = 0; i < np - 1; ++i)
2561 {
2562 for (int j = 0; j < np - 1; ++j)
2563 {
2564 rowp1 += np;
2565 for (int k = 0; k < np - 1; ++k)
2566 {
2567 conn[cnt++] = plane + row + k;
2568 conn[cnt++] = plane + row + k + 1;
2569 conn[cnt++] = plane + rowp1 + k;
2570
2571 conn[cnt++] = plane + rowp1 + k + 1;
2572 conn[cnt++] = plane + rowp1 + k;
2573 conn[cnt++] = plane + row + k + 1;
2574 }
2575 row += np;
2576 }
2577 plane += np * np;
2578 }
2579}

References Nektar::StdRegions::StdExpansion::m_base.

◆ v_GetTraceBasisKey()

const LibUtilities::BasisKey Nektar::StdRegions::StdHexExp::v_GetTraceBasisKey ( const int  i,
const int  k 
) const
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 745 of file StdHexExp.cpp.

747{
748 ASSERTL2(i >= 0 && i <= 5, "face id is out of range");
749 ASSERTL2(k >= 0 && k <= 1, "basis key id is out of range");
750
751 int dir = k;
752 switch (i)
753 {
754 case 0:
755 case 5:
756 dir = k;
757 break;
758 case 1:
759 case 3:
760 dir = 2 * k;
761 break;
762 case 2:
763 case 4:
764 dir = k + 1;
765 break;
766 }
767
769 m_base[dir]->GetNumPoints(),
770 m_base[dir]->GetNumModes());
771}
LibUtilities::BasisKey EvaluateQuadFaceBasisKey(const int facedir, const LibUtilities::BasisType faceDirBasisType, const int numpoints, const int nummodes)

References ASSERTL2, Nektar::StdRegions::EvaluateQuadFaceBasisKey(), Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::GetNumPoints(), 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 1137 of file StdHexExp.cpp.

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

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 1849 of file StdHexExp.cpp.

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

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, CellMLToNektar.cellml_metadata::p, Nektar::UnitTests::q(), and v_GetTraceIntNcoeffs().

◆ v_GetTraceIntNcoeffs()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 680 of file StdHexExp.cpp.

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

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 663 of file StdHexExp.cpp.

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

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 801 of file StdHexExp.cpp.

803{
804 int nummodes[3] = {m_base[0]->GetNumModes(), m_base[1]->GetNumModes(),
805 m_base[2]->GetNumModes()};
806 switch (fid)
807 {
808 case 0:
809 case 5:
810 {
811 numModes0 = nummodes[0];
812 numModes1 = nummodes[1];
813 }
814 break;
815 case 1:
816 case 3:
817 {
818 numModes0 = nummodes[0];
819 numModes1 = nummodes[2];
820 }
821 break;
822 case 2:
823 case 4:
824 {
825 numModes0 = nummodes[1];
826 numModes1 = nummodes[2];
827 }
828 break;
829 default:
830 {
831 ASSERTL0(false, "fid out of range");
832 }
833 break;
834 }
835
836 if (faceOrient >= eDir1FwdDir2_Dir2FwdDir1)
837 {
838 std::swap(numModes0, numModes1);
839 }
840}

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 697 of file StdHexExp.cpp.

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

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 715 of file StdHexExp.cpp.

717{
718 ASSERTL2(i >= 0 && i <= 5, "face id is out of range");
719 ASSERTL2(j == 0 || j == 1, "face direction is out of range");
720
721 if (i == 0 || i == 5)
722 {
723 return m_base[j]->GetPointsKey();
724 }
725 else if (i == 1 || i == 3)
726 {
727 return m_base[2 * j]->GetPointsKey();
728 }
729 else
730 {
731 return m_base[j + 1]->GetPointsKey();
732 }
733}

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 849 of file StdHexExp.cpp.

850{
853 "BasisType is not a boundary interior form");
856 "BasisType is not a boundary interior form");
859 "BasisType is not a boundary interior form");
860
861 ASSERTL1((localVertexId >= 0) && (localVertexId < 8),
862 "local vertex id must be between 0 and 7");
863
864 int p = 0;
865 int q = 0;
866 int r = 0;
867
868 // Retrieve the number of modes in each dimension.
869 int nummodes[3] = {m_base[0]->GetNumModes(), m_base[1]->GetNumModes(),
870 m_base[2]->GetNumModes()};
871
872 if (useCoeffPacking == true) // follow packing of coefficients i.e q,r,p
873 {
874 if (localVertexId > 3)
875 {
877 {
878 r = nummodes[2] - 1;
879 }
880 else
881 {
882 r = 1;
883 }
884 }
885
886 switch (localVertexId % 4)
887 {
888 case 0:
889 break;
890 case 1:
891 {
893 {
894 p = nummodes[0] - 1;
895 }
896 else
897 {
898 p = 1;
899 }
900 }
901 break;
902 case 2:
903 {
905 {
906 q = nummodes[1] - 1;
907 }
908 else
909 {
910 q = 1;
911 }
912 }
913 break;
914 case 3:
915 {
917 {
918 p = nummodes[0] - 1;
919 q = nummodes[1] - 1;
920 }
921 else
922 {
923 p = 1;
924 q = 1;
925 }
926 }
927 break;
928 }
929 }
930 else
931 {
932 // Right face (vertices 1,2,5,6)
933 if ((localVertexId % 4) % 3 > 0)
934 {
936 {
937 p = nummodes[0] - 1;
938 }
939 else
940 {
941 p = 1;
942 }
943 }
944 // Back face (vertices 2,3,6,7)
945 if (localVertexId % 4 > 1)
946 {
948 {
949 q = nummodes[1] - 1;
950 }
951 else
952 {
953 q = 1;
954 }
955 }
956
957 // Top face (vertices 4,5,6,7)
958 if (localVertexId > 3)
959 {
961 {
962 r = nummodes[2] - 1;
963 }
964 else
965 {
966 r = 1;
967 }
968 }
969 }
970 // Compute the local number.
971 return r * nummodes[0] * nummodes[1] + q * nummodes[0] + p;
972}

References ASSERTL1, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::m_base, CellMLToNektar.cellml_metadata::p, and Nektar::UnitTests::q().

◆ 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.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 2314 of file StdHexExp.cpp.

2317{
2318 StdHexExp::v_HelmholtzMatrixOp_MatFree(inarray, outarray, mkey);
2319}
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override

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

◆ v_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.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 334 of file StdHexExp.cpp.

336{
337 if (m_base[0]->Collocation() && m_base[1]->Collocation() &&
338 m_base[2]->Collocation())
339 {
340 MultiplyByQuadratureMetric(inarray, outarray);
341 }
342 else
343 {
344 StdHexExp::v_IProductWRTBase_SumFac(inarray, outarray);
345 }
346}
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:729
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multbyweights=true) override
Definition: StdHexExp.cpp:351

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.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 351 of file StdHexExp.cpp.

354{
355 int nquad0 = m_base[0]->GetNumPoints();
356 int nquad1 = m_base[1]->GetNumPoints();
357 int nquad2 = m_base[2]->GetNumPoints();
358 int order0 = m_base[0]->GetNumModes();
359 int order1 = m_base[1]->GetNumModes();
360
361 Array<OneD, NekDouble> wsp(nquad0 * nquad1 * (nquad2 + order0) +
362 order0 * order1 * nquad2);
363
364 if (multiplybyweights)
365 {
366 Array<OneD, NekDouble> tmp(inarray.size());
367 MultiplyByQuadratureMetric(inarray, tmp);
368
370 m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
371 tmp, outarray, wsp, true, true, true);
372 }
373 else
374 {
376 m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
377 inarray, outarray, wsp, true, true, true);
378 }
379}
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 385 of file StdHexExp.cpp.

392{
393 int nquad0 = m_base[0]->GetNumPoints();
394 int nquad1 = m_base[1]->GetNumPoints();
395 int nquad2 = m_base[2]->GetNumPoints();
396 int nmodes0 = m_base[0]->GetNumModes();
397 int nmodes1 = m_base[1]->GetNumModes();
398 int nmodes2 = m_base[2]->GetNumModes();
399
400 bool colldir0 = doCheckCollDir0 ? (m_base[0]->Collocation()) : false;
401 bool colldir1 = doCheckCollDir1 ? (m_base[1]->Collocation()) : false;
402 bool colldir2 = doCheckCollDir2 ? (m_base[2]->Collocation()) : false;
403
404 if (colldir0 && colldir1 && colldir2)
405 {
406 Vmath::Vcopy(m_ncoeffs, inarray.get(), 1, outarray.get(), 1);
407 }
408 else
409 {
410 ASSERTL1(wsp.size() >= nmodes0 * nquad2 * (nquad1 + nmodes1),
411 "Insufficient workspace size");
412
413 Array<OneD, NekDouble> tmp0 = wsp;
414 Array<OneD, NekDouble> tmp1 = wsp + nmodes0 * nquad1 * nquad2;
415
416 if (colldir0)
417 {
418 // reshuffle data for next operation.
419 for (int n = 0; n < nmodes0; ++n)
420 {
421 Vmath::Vcopy(nquad1 * nquad2, inarray.get() + n, nquad0,
422 tmp0.get() + nquad1 * nquad2 * n, 1);
423 }
424 }
425 else
426 {
427 Blas::Dgemm('T', 'N', nquad1 * nquad2, nmodes0, nquad0, 1.0,
428 inarray.get(), nquad0, base0.get(), nquad0, 0.0,
429 tmp0.get(), nquad1 * nquad2);
430 }
431
432 if (colldir1)
433 {
434 // reshuffle data for next operation.
435 for (int n = 0; n < nmodes1; ++n)
436 {
437 Vmath::Vcopy(nquad2 * nmodes0, tmp0.get() + n, nquad1,
438 tmp1.get() + nquad2 * nmodes0 * n, 1);
439 }
440 }
441 else
442 {
443 Blas::Dgemm('T', 'N', nquad2 * nmodes0, nmodes1, nquad1, 1.0,
444 tmp0.get(), nquad1, base1.get(), nquad1, 0.0,
445 tmp1.get(), nquad2 * nmodes0);
446 }
447
448 if (colldir2)
449 {
450 // reshuffle data for next operation.
451 for (int n = 0; n < nmodes2; ++n)
452 {
453 Vmath::Vcopy(nmodes0 * nmodes1, tmp1.get() + n, nquad2,
454 outarray.get() + nmodes0 * nmodes1 * n, 1);
455 }
456 }
457 else
458 {
459 Blas::Dgemm('T', 'N', nmodes0 * nmodes1, nmodes2, nquad2, 1.0,
460 tmp1.get(), nquad2, base2.get(), nquad2, 0.0,
461 outarray.get(), nmodes0 * nmodes1);
462 }
463 }
464}

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.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 466 of file StdHexExp.cpp.

469{
470 StdHexExp::IProductWRTDerivBase_SumFac(dir, inarray, outarray);
471}
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.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 473 of file StdHexExp.cpp.

476{
477 ASSERTL0((dir == 0) || (dir == 1) || (dir == 2),
478 "input dir is out of range");
479
480 int nquad1 = m_base[1]->GetNumPoints();
481 int nquad2 = m_base[2]->GetNumPoints();
482 int order0 = m_base[0]->GetNumModes();
483 int order1 = m_base[1]->GetNumModes();
484
485 // If outarray > inarray then no need for temporary storage.
486 Array<OneD, NekDouble> tmp = outarray;
487 if (outarray.size() < inarray.size())
488 {
489 tmp = Array<OneD, NekDouble>(inarray.size());
490 }
491
492 // Need workspace for sumfackernel though
493 Array<OneD, NekDouble> wsp(order0 * nquad2 * (nquad1 + order1));
494
495 // multiply by integration constants
496 MultiplyByQuadratureMetric(inarray, tmp);
497
498 // perform sum-factorisation
499 switch (dir)
500 {
501 case 0:
503 m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
504 m_base[2]->GetBdata(), tmp, outarray, wsp, false, true, true);
505 break;
506 case 1:
508 m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
509 m_base[2]->GetBdata(), tmp, outarray, wsp, true, false, true);
510 break;
511 case 2:
513 m_base[0]->GetBdata(), m_base[1]->GetBdata(),
514 m_base[2]->GetDbdata(), tmp, outarray, wsp, true, true, false);
515 break;
516 }
517}

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.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 2291 of file StdHexExp.cpp.

2294{
2295 StdHexExp::v_LaplacianMatrixOp_MatFree(inarray, outarray, mkey);
2296}
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey) override

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

◆ v_LaplacianMatrixOp() [2/2]

void Nektar::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.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 2298 of file StdHexExp.cpp.

2302{
2303 StdExpansion::LaplacianMatrixOp_MatFree(k1, k2, inarray, outarray, mkey);
2304}
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 527 of file StdHexExp.cpp.

529{
530 xi[0] = eta[0];
531 xi[1] = eta[1];
532 xi[2] = eta[2];
533}

◆ 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 519 of file StdHexExp.cpp.

521{
522 eta[0] = xi[0];
523 eta[1] = xi[1];
524 eta[2] = xi[2];
525}

◆ 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.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 2284 of file StdHexExp.cpp.

2287{
2288 StdExpansion::MassMatrixOp_MatFree(inarray, outarray, mkey);
2289}
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 2321 of file StdHexExp.cpp.

2324{
2325 int nquad0 = m_base[0]->GetNumPoints();
2326 int nquad1 = m_base[1]->GetNumPoints();
2327 int nquad2 = m_base[2]->GetNumPoints();
2328 int nq01 = nquad0 * nquad1;
2329 int nq12 = nquad1 * nquad2;
2330
2331 const Array<OneD, const NekDouble> &w0 = m_base[0]->GetW();
2332 const Array<OneD, const NekDouble> &w1 = m_base[1]->GetW();
2333 const Array<OneD, const NekDouble> &w2 = m_base[2]->GetW();
2334
2335 for (int i = 0; i < nq12; ++i)
2336 {
2337 Vmath::Vmul(nquad0, inarray.get() + i * nquad0, 1, w0.get(), 1,
2338 outarray.get() + i * nquad0, 1);
2339 }
2340
2341 for (int i = 0; i < nq12; ++i)
2342 {
2343 Vmath::Smul(nquad0, w1[i % nquad1], outarray.get() + i * nquad0, 1,
2344 outarray.get() + i * nquad0, 1);
2345 }
2346
2347 for (int i = 0; i < nquad2; ++i)
2348 {
2349 Vmath::Smul(nq01, w2[i], outarray.get() + i * nq01, 1,
2350 outarray.get() + i * nq01, 1);
2351 }
2352}

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 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 4 * (nmodes0 + nmodes1 + nmodes2) + 8);
642}

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 644 of file StdHexExp.cpp.

645{
648 "BasisType is not a boundary interior form");
651 "BasisType is not a boundary interior form");
654 "BasisType is not a boundary interior form");
655
656 int nmodes0 = m_base[0]->GetNumModes();
657 int nmodes1 = m_base[1]->GetNumModes();
658 int nmodes2 = m_base[2]->GetNumModes();
659
660 return 2 * (nmodes0 * nmodes1 + nmodes0 * nmodes2 + nmodes1 * nmodes2);
661}

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.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 86 of file StdHexExp.cpp.

90{
91 StdExpansion3D::PhysTensorDeriv(inarray, out_d0, out_d1, out_d2);
92}
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.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 100 of file StdHexExp.cpp.

103{
104 switch (dir)
105 {
106 case 0:
107 {
108 PhysDeriv(inarray, outarray, NullNekDouble1DArray,
110 }
111 break;
112 case 1:
113 {
114 PhysDeriv(inarray, NullNekDouble1DArray, outarray,
116 }
117 break;
118 case 2:
119 {
121 outarray);
122 }
123 break;
124 default:
125 {
126 ASSERTL1(false, "input dir is out of range");
127 }
128 break;
129 }
130}
void PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Definition: StdExpansion.h:855
static Array< OneD, NekDouble > NullNekDouble1DArray

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

◆ v_PhysEvaluate()

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

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 1127 of file StdHexExp.cpp.

1130{
1131 return BaryTensorDeriv(coord, inarray, firstOrderDerivs);
1132}
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 
)
finaloverrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 583 of file StdHexExp.cpp.

585{
586 ASSERTL2(coords[0] > -1 - NekConstants::kNekZeroTol, "coord[0] < -1");
587 ASSERTL2(coords[0] < 1 + NekConstants::kNekZeroTol, "coord[0] > 1");
588 ASSERTL2(coords[1] > -1 - NekConstants::kNekZeroTol, "coord[1] < -1");
589 ASSERTL2(coords[1] < 1 + NekConstants::kNekZeroTol, "coord[1] > 1");
590 ASSERTL2(coords[2] > -1 - NekConstants::kNekZeroTol, "coord[2] < -1");
591 ASSERTL2(coords[2] < 1 + NekConstants::kNekZeroTol, "coord[2] > 1");
592
593 const int nm0 = m_base[0]->GetNumModes();
594 const int nm1 = m_base[1]->GetNumModes();
595 const int mode2 = mode / (nm0 * nm1);
596 const int mode1 = (mode - mode2 * nm0 * nm1) / nm0;
597 const int mode0 = (mode - mode2 * nm0 * nm1) % nm0;
598
599 return StdExpansion::BaryEvaluateBasis<0>(coords[0], mode0) *
600 StdExpansion::BaryEvaluateBasis<1>(coords[1], mode1) *
601 StdExpansion::BaryEvaluateBasis<2>(coords[2], mode2);
602}
static const NekDouble kNekZeroTol

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

◆ v_StdPhysDeriv() [1/2]

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 132 of file StdHexExp.cpp.

136{
137 StdHexExp::v_PhysDeriv(inarray, out_d0, out_d1, out_d2);
138}
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:86

References v_PhysDeriv().

◆ v_StdPhysDeriv() [2/2]

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 140 of file StdHexExp.cpp.

143{
144 StdHexExp::v_PhysDeriv(dir, inarray, outarray);
145}

References v_PhysDeriv().

◆ v_SVVLaplacianFilter()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 2354 of file StdHexExp.cpp.

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

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, and Nektar::StdRegions::StdExpansion::m_base.

◆ 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.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 2306 of file StdHexExp.cpp.

2310{
2311 StdExpansion::WeakDerivMatrixOp_MatFree(i, inarray, outarray, mkey);
2312}
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().