Nektar++
Public Member Functions | Protected Member Functions | Private Member Functions | List of all members
Nektar::StdRegions::StdTetExp Class Reference

#include <StdTetExp.h>

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

Public Member Functions

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

Protected Member Functions

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

Private Member Functions

int GetMode (const int i, const int j, const int k)
 Compute the mode number in the expansion for a particular tensorial combination. More...
 

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

Definition at line 43 of file StdTetExp.h.

Constructor & Destructor Documentation

◆ StdTetExp() [1/3]

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

Definition at line 42 of file StdTetExp.cpp.

46 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
47 3, Ba, Bb, Bc),
49 Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
50 Ba, Bb, Bc)
51{
52 ASSERTL0(Ba.GetNumModes() <= Bb.GetNumModes(),
53 "order in 'a' direction is higher than order "
54 "in 'b' direction");
55 ASSERTL0(Ba.GetNumModes() <= Bc.GetNumModes(),
56 "order in 'a' direction is higher than order "
57 "in 'c' direction");
58 ASSERTL0(Bb.GetNumModes() <= Bc.GetNumModes(),
59 "order in 'b' direction is higher than order "
60 "in 'c' direction");
61}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
StdExpansion()
Default Constructor.
int getNumberOfCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:187

References ASSERTL0, and Nektar::LibUtilities::BasisKey::GetNumModes().

◆ StdTetExp() [2/3]

Nektar::StdRegions::StdTetExp::StdTetExp ( const LibUtilities::BasisKey Ba,
const LibUtilities::BasisKey Bb,
const LibUtilities::BasisKey Bc,
NekDouble coeffs,
NekDouble phys 
)

◆ StdTetExp() [3/3]

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

◆ ~StdTetExp()

Nektar::StdRegions::StdTetExp::~StdTetExp ( )
overridedefault

Member Function Documentation

◆ DetShapeType()

LibUtilities::ShapeType Nektar::StdRegions::StdTetExp::DetShapeType ( ) const
inline

◆ GetMode()

int Nektar::StdRegions::StdTetExp::GetMode ( const int  I,
const int  J,
const int  K 
)
private

Compute the mode number in the expansion for a particular tensorial combination.

Modes are numbered with the r index travelling fastest, followed by q and then p, and each q-r plane is of size (Q+1)*(Q+2)/2+max(0,R-Q-p)*Q. For example, when P=2, Q=3 and R=4 (nm0=3, nm1 = 4, nm2 = 5) the indexing inside each q-r plane (with r increasing upwards and q to the right) is:

4 3 8 17 2 7 11 16 20 25 1 6 10 13 15 19 22 24 27 0 5 9 12 14 18 21 23 26

Geometrically they can be interpreted as

p = 0: p = 2: p = 1:

1 4 8 17 3 11 7 25 16 20 2 10 13 6 24 27 15 19 22 0 9 12 5 23 26 14 18 21

so we have the following breakdown

Vertices V[0,1,2,3] = [0, 14, 5, 1] Edges E[0,1,2,3,4,5,6] =[[23],[18, 21], [9, 12], [2,3,4], [15, 16, 17], [6,7,8]] Faces F[0.1,2,3] = [[26], [24,25], [19, 22, 20], [10, 13, 11] Interior [27] Note that in this element, we must have that \( P \leq Q \leq R\).

Definition at line 1918 of file StdTetExp.cpp.

1919{
1920 const int Q = m_base[1]->GetNumModes();
1921 const int R = m_base[2]->GetNumModes();
1922
1923 int i, j, q_hat, k_hat;
1924 int cnt = 0;
1925
1926 // Traverse to q-r plane number I
1927 for (i = 0; i < I; ++i)
1928 {
1929 // Size of triangle part
1930 q_hat = Q - i;
1931 // Size of rectangle part
1932 k_hat = R - Q;
1933 cnt += q_hat * (q_hat + 1) / 2 + k_hat * (Q - i);
1934 }
1935
1936 // Traverse to q column J
1937 q_hat = R - I;
1938 for (j = 0; j < J; ++j)
1939 {
1940 cnt += q_hat;
1941 q_hat--;
1942 }
1943
1944 // Traverse up stacks to K
1945 cnt += K;
1946
1947 return cnt;
1948}
Array< OneD, LibUtilities::BasisSharedPtr > m_base

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

Referenced by v_GetBoundaryMap(), v_GetEdgeInteriorToElementMap(), v_GetInteriorMap(), v_GetTraceCoeffMap(), v_GetTraceInteriorToElementMap(), and v_GetVertexMap().

◆ PhysEvaluate3D()

NekDouble Nektar::StdRegions::StdTetExp::PhysEvaluate3D ( const Array< OneD, const NekDouble > &  coords,
const Array< OneD, const NekDouble > &  physvals 
)

Single Point Evaluation.

◆ v_BwdTrans()

void Nektar::StdRegions::StdTetExp::v_BwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual
Note
'r' (base[2]) runs fastest in this element

\( u^{\delta} (\xi_{1i}, \xi_{2j}, \xi_{3k}) = \sum_{m(pqr)} \hat u_{pqr} \phi_{pqr} (\xi_{1i}, \xi_{2j}, \xi_{3k})\)

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_{pq}^b (\xi_{2j}) \lbrace { \sum_{r=0}^{Q_z} \hat u_{pqr} \psi_{pqr}^c (\xi_{3k}) \rbrace} \rbrace}. \) And sumfactorizing step of the form is as:\

\( f_{pq} (\xi_{3k}) = \sum_{r=0}^{Q_z} \hat u_{pqr} \psi_{pqr}^c (\xi_{3k}),\\ g_{p} (\xi_{2j}, \xi_{3k}) = \sum_{r=0}^{Q_y} \psi_{pq}^b (\xi_{2j}) f_{pq} (\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}). \)

Implements Nektar::StdRegions::StdExpansion.

Definition at line 266 of file StdTetExp.cpp.

268{
271 "Basis[1] is not a general tensor type");
272
275 "Basis[2] is not a general tensor type");
276
277 if (m_base[0]->Collocation() && m_base[1]->Collocation() &&
278 m_base[2]->Collocation())
279 {
281 m_base[2]->GetNumPoints(),
282 inarray, 1, outarray, 1);
283 }
284 else
285 {
286 StdTetExp::v_BwdTrans_SumFac(inarray, outarray);
287 }
288}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:156
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:218
void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdTetExp.cpp:293
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:49
@ eModified_C
Principle Modified Functions .
Definition: BasisType.h:50
@ eOrtho_C
Principle Orthogonal Functions .
Definition: BasisType.h:46
@ eOrtho_B
Principle Orthogonal Functions .
Definition: BasisType.h:44
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References ASSERTL1, Nektar::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, v_BwdTrans_SumFac(), and Vmath::Vcopy().

Referenced by v_FillMode().

◆ v_BwdTrans_SumFac()

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

Sum-factorisation implementation of the BwdTrans operation.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 293 of file StdTetExp.cpp.

295{
296 int nquad1 = m_base[1]->GetNumPoints();
297 int nquad2 = m_base[2]->GetNumPoints();
298 int order0 = m_base[0]->GetNumModes();
299 int order1 = m_base[1]->GetNumModes();
300
301 Array<OneD, NekDouble> wsp(nquad2 * order0 * (2 * order1 - order0 + 1) / 2 +
302 nquad2 * nquad1 * order0);
303
304 BwdTrans_SumFacKernel(m_base[0]->GetBdata(), m_base[1]->GetBdata(),
305 m_base[2]->GetBdata(), inarray, outarray, wsp, true,
306 true, true);
307}
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(), and Nektar::StdRegions::StdExpansion::m_base.

Referenced by v_BwdTrans(), and Nektar::StdRegions::StdNodalTetExp::v_BwdTrans_SumFac().

◆ v_BwdTrans_SumFacKernel()

void Nektar::StdRegions::StdTetExp::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 322 of file StdTetExp.cpp.

331{
332 int nquad0 = m_base[0]->GetNumPoints();
333 int nquad1 = m_base[1]->GetNumPoints();
334 int nquad2 = m_base[2]->GetNumPoints();
335
336 int order0 = m_base[0]->GetNumModes();
337 int order1 = m_base[1]->GetNumModes();
338 int order2 = m_base[2]->GetNumModes();
339
340 Array<OneD, NekDouble> tmp = wsp;
341 Array<OneD, NekDouble> tmp1 =
342 tmp + nquad2 * order0 * (2 * order1 - order0 + 1) / 2;
343
344 int i, j, mode, mode1, cnt;
345
346 // Perform summation over '2' direction
347 mode = mode1 = cnt = 0;
348 for (i = 0; i < order0; ++i)
349 {
350 for (j = 0; j < order1 - i; ++j, ++cnt)
351 {
352 Blas::Dgemv('N', nquad2, order2 - i - j, 1.0,
353 base2.get() + mode * nquad2, nquad2,
354 inarray.get() + mode1, 1, 0.0, tmp.get() + cnt * nquad2,
355 1);
356 mode += order2 - i - j;
357 mode1 += order2 - i - j;
358 }
359 // increment mode in case order1!=order2
360 for (j = order1 - i; j < order2 - i; ++j)
361 {
362 mode += order2 - i - j;
363 }
364 }
365
366 // fix for modified basis by adding split of top singular
367 // vertex mode - currently (1+c)/2 x (1-b)/2 x (1-a)/2
368 // component is evaluated
370 {
371 // top singular vertex - (1+c)/2 x (1+b)/2 x (1-a)/2 component
372 Blas::Daxpy(nquad2, inarray[1], base2.get() + nquad2, 1,
373 &tmp[0] + nquad2, 1);
374
375 // top singular vertex - (1+c)/2 x (1-b)/2 x (1+a)/2 component
376 Blas::Daxpy(nquad2, inarray[1], base2.get() + nquad2, 1,
377 &tmp[0] + order1 * nquad2, 1);
378 }
379
380 // Perform summation over '1' direction
381 mode = 0;
382 for (i = 0; i < order0; ++i)
383 {
384 Blas::Dgemm('N', 'T', nquad1, nquad2, order1 - i, 1.0,
385 base1.get() + mode * nquad1, nquad1,
386 tmp.get() + mode * nquad2, nquad2, 0.0,
387 tmp1.get() + i * nquad1 * nquad2, nquad1);
388 mode += order1 - i;
389 }
390
391 // fix for modified basis by adding additional split of
392 // top and base singular vertex modes as well as singular
393 // edge
395 {
396 // use tmp to sort out singular vertices and
397 // singular edge components with (1+b)/2 (1+a)/2 form
398 for (i = 0; i < nquad2; ++i)
399 {
400 Blas::Daxpy(nquad1, tmp[nquad2 + i], base1.get() + nquad1, 1,
401 &tmp1[nquad1 * nquad2] + i * nquad1, 1);
402 }
403 }
404
405 // Perform summation over '0' direction
406 Blas::Dgemm('N', 'T', nquad0, nquad1 * nquad2, order0, 1.0, base0.get(),
407 nquad0, tmp1.get(), nquad1 * nquad2, 0.0, outarray.get(),
408 nquad0);
409}
static void Dgemv(const char &trans, const int &m, const int &n, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = alpha A x plus beta y where A[m x n].
Definition: Blas.hpp:211
static void Dgemm(const char &transa, const char &transb, const int &m, const int &n, const int &k, const double &alpha, const double *a, const int &lda, const double *b, const int &ldb, const double &beta, double *c, const int &ldc)
BLAS level 3: Matrix-matrix multiply C = A x B where op(A)[m x k], op(B)[k x n], C[m x n] DGEMM perfo...
Definition: Blas.hpp:383
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
Definition: Blas.hpp:135
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48

References Blas::Daxpy(), Blas::Dgemm(), Blas::Dgemv(), Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), and Nektar::StdRegions::StdExpansion::m_base.

◆ v_CalcNumberOfCoefficients()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1121 of file StdTetExp.cpp.

1123{
1125 nummodes[modes_offset], nummodes[modes_offset + 1],
1126 nummodes[modes_offset + 2]);
1127 modes_offset += 3;
1128
1129 return nmodes;
1130}

References Nektar::LibUtilities::StdTetData::getNumberOfCoefficients().

◆ v_CreateStdMatrix()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TetExp.

Definition at line 1874 of file StdTetExp.cpp.

1875{
1876 return v_GenMatrix(mkey);
1877}
DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override
Definition: StdTetExp.cpp:1786

References v_GenMatrix().

◆ v_DetShapeType()

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

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TetExp.

Definition at line 970 of file StdTetExp.cpp.

971{
973}

References Nektar::LibUtilities::eTetrahedron.

◆ v_FillMode()

void Nektar::StdRegions::StdTetExp::v_FillMode ( const int  mode,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 808 of file StdTetExp.cpp.

809{
810 Array<OneD, NekDouble> tmp(m_ncoeffs, 0.0);
811 tmp[mode] = 1.0;
812 StdTetExp::v_BwdTrans(tmp, outarray);
813}
void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdTetExp.cpp:266

References Nektar::StdRegions::StdExpansion::m_ncoeffs, and v_BwdTrans().

Referenced by Nektar::StdRegions::StdNodalTetExp::GenNBasisTransMatrix().

◆ v_FwdTrans()

void Nektar::StdRegions::StdTetExp::v_FwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual
Parameters
inarrayarray of physical quadrature points to be transformed.
outarrayupdated array of expansion coefficients.

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TetExp.

Definition at line 416 of file StdTetExp.cpp.

418{ // int numMax = nmodes0;
419 v_IProductWRTBase(inarray, outarray);
420
421 // get Mass matrix inverse
422 StdMatrixKey masskey(eInvMass, DetShapeType(), *this);
423 DNekMatSharedPtr matsys = GetStdMatrix(masskey);
424
425 // copy inarray in case inarray == outarray
426 DNekVec in(m_ncoeffs, outarray);
427 DNekVec out(m_ncoeffs, outarray, eWrapper);
428
429 out = (*matsys) * in;
430}
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:603
void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdTetExp.cpp:466
LibUtilities::ShapeType DetShapeType() const
Definition: StdTetExp.h:56
NekVector< NekDouble > DNekVec
Definition: NekTypeDefs.hpp:48
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75

References DetShapeType(), Nektar::StdRegions::eInvMass, Nektar::eWrapper, Nektar::StdRegions::StdExpansion::GetStdMatrix(), Nektar::StdRegions::StdExpansion::m_ncoeffs, and v_IProductWRTBase().

◆ v_GenMatrix()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TetExp.

Definition at line 1786 of file StdTetExp.cpp.

1787{
1788
1789 MatrixType mtype = mkey.GetMatrixType();
1790
1791 DNekMatSharedPtr Mat;
1792
1793 switch (mtype)
1794 {
1796 {
1797 int nq0 = m_base[0]->GetNumPoints();
1798 int nq1 = m_base[1]->GetNumPoints();
1799 int nq2 = m_base[2]->GetNumPoints();
1800 int nq;
1801
1802 // take definition from key
1803 if (mkey.ConstFactorExists(eFactorConst))
1804 {
1805 nq = (int)mkey.GetConstFactor(eFactorConst);
1806 }
1807 else
1808 {
1809 nq = max(nq0, max(nq1, nq2));
1810 }
1811
1812 int neq =
1814 Array<OneD, Array<OneD, NekDouble>> coords(neq);
1815 Array<OneD, NekDouble> coll(3);
1816 Array<OneD, DNekMatSharedPtr> I(3);
1817 Array<OneD, NekDouble> tmp(nq0);
1818
1819 Mat =
1820 MemoryManager<DNekMat>::AllocateSharedPtr(neq, nq0 * nq1 * nq2);
1821 int cnt = 0;
1822
1823 for (int i = 0; i < nq; ++i)
1824 {
1825 for (int j = 0; j < nq - i; ++j)
1826 {
1827 for (int k = 0; k < nq - i - j; ++k, ++cnt)
1828 {
1829 coords[cnt] = Array<OneD, NekDouble>(3);
1830 coords[cnt][0] = -1.0 + 2 * k / (NekDouble)(nq - 1);
1831 coords[cnt][1] = -1.0 + 2 * j / (NekDouble)(nq - 1);
1832 coords[cnt][2] = -1.0 + 2 * i / (NekDouble)(nq - 1);
1833 }
1834 }
1835 }
1836
1837 for (int i = 0; i < neq; ++i)
1838 {
1839 LocCoordToLocCollapsed(coords[i], coll);
1840
1841 I[0] = m_base[0]->GetI(coll);
1842 I[1] = m_base[1]->GetI(coll + 1);
1843 I[2] = m_base[2]->GetI(coll + 2);
1844
1845 // interpolate first coordinate direction
1846 NekDouble fac;
1847 for (int k = 0; k < nq2; ++k)
1848 {
1849 for (int j = 0; j < nq1; ++j)
1850 {
1851
1852 fac = (I[1]->GetPtr())[j] * (I[2]->GetPtr())[k];
1853 Vmath::Smul(nq0, fac, I[0]->GetPtr(), 1, tmp, 1);
1854
1855 Vmath::Vcopy(nq0, &tmp[0], 1,
1856 Mat->GetRawPtr() + k * nq0 * nq1 * neq +
1857 j * nq0 * neq + i,
1858 neq);
1859 }
1860 }
1861 }
1862 }
1863 break;
1864 default:
1865 {
1867 }
1868 break;
1869 }
1870
1871 return Mat;
1872}
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
double NekDouble
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100

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

Referenced by v_CreateStdMatrix().

◆ v_GetBoundaryMap()

void Nektar::StdRegions::StdTetExp::v_GetBoundaryMap ( Array< OneD, unsigned int > &  outarray)
overrideprotectedvirtual

List of all boundary modes in the the expansion.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1318 of file StdTetExp.cpp.

1319{
1322 "BasisType is not a boundary interior form");
1325 "BasisType is not a boundary interior form");
1328 "BasisType is not a boundary interior form");
1329
1330 int P = m_base[0]->GetNumModes();
1331 int Q = m_base[1]->GetNumModes();
1332 int R = m_base[2]->GetNumModes();
1333
1334 int i, j, k;
1335 int idx = 0;
1336
1337 int nBnd = NumBndryCoeffs();
1338
1339 if (outarray.size() != nBnd)
1340 {
1341 outarray = Array<OneD, unsigned int>(nBnd);
1342 }
1343
1344 for (i = 0; i < P; ++i)
1345 {
1346 // First two Q-R planes are entirely boundary modes
1347 if (i < 2)
1348 {
1349 for (j = 0; j < Q - i; j++)
1350 {
1351 for (k = 0; k < R - i - j; ++k)
1352 {
1353 outarray[idx++] = GetMode(i, j, k);
1354 }
1355 }
1356 }
1357 // Remaining Q-R planes contain boundary modes on bottom and
1358 // left edge.
1359 else
1360 {
1361 for (k = 0; k < R - i; ++k)
1362 {
1363 outarray[idx++] = GetMode(i, 0, k);
1364 }
1365 for (j = 1; j < Q - i; ++j)
1366 {
1367 outarray[idx++] = GetMode(i, j, 0);
1368 }
1369 }
1370 }
1371}
int GetMode(const int i, const int j, const int k)
Compute the mode number in the expansion for a particular tensorial combination.
Definition: StdTetExp.cpp:1918
@ P
Monomial polynomials .
Definition: BasisType.h:62
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:56

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

◆ v_GetCoords()

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

Definition at line 1158 of file StdTetExp.cpp.

1161{
1162 Array<OneD, const NekDouble> eta_x = m_base[0]->GetZ();
1163 Array<OneD, const NekDouble> eta_y = m_base[1]->GetZ();
1164 Array<OneD, const NekDouble> eta_z = m_base[2]->GetZ();
1165 int Qx = GetNumPoints(0);
1166 int Qy = GetNumPoints(1);
1167 int Qz = GetNumPoints(2);
1168
1169 // Convert collapsed coordinates into cartesian coordinates: eta
1170 // --> xi
1171 for (int k = 0; k < Qz; ++k)
1172 {
1173 for (int j = 0; j < Qy; ++j)
1174 {
1175 for (int i = 0; i < Qx; ++i)
1176 {
1177 int s = i + Qx * (j + Qy * k);
1178 xi_x[s] =
1179 (eta_x[i] + 1.0) * (1.0 - eta_y[j]) * (1.0 - eta_z[k]) / 4 -
1180 1.0;
1181 xi_y[s] = (eta_y[j] + 1.0) * (1.0 - eta_z[k]) / 2 - 1.0;
1182 xi_z[s] = eta_z[k];
1183 }
1184 }
1185 }
1186}

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

◆ v_GetEdgeInteriorToElementMap()

void Nektar::StdRegions::StdTetExp::v_GetEdgeInteriorToElementMap ( const int  eid,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray,
const Orientation  edgeOrient = eDir1FwdDir1_Dir2FwdDir2 
)
overrideprotectedvirtual

Maps interior modes of an edge to the elemental modes.

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 1580 of file StdTetExp.cpp.

1583{
1584 int i;
1585 const int P = m_base[0]->GetNumModes();
1586 const int Q = m_base[1]->GetNumModes();
1587 const int R = m_base[2]->GetNumModes();
1588
1589 const int nEdgeIntCoeffs = v_GetEdgeNcoeffs(eid) - 2;
1590
1591 if (maparray.size() != nEdgeIntCoeffs)
1592 {
1593 maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
1594 }
1595 else
1596 {
1597 fill(maparray.get(), maparray.get() + nEdgeIntCoeffs, 0);
1598 }
1599
1600 if (signarray.size() != nEdgeIntCoeffs)
1601 {
1602 signarray = Array<OneD, int>(nEdgeIntCoeffs, 1);
1603 }
1604 else
1605 {
1606 fill(signarray.get(), signarray.get() + nEdgeIntCoeffs, 1);
1607 }
1608
1609 switch (eid)
1610 {
1611 case 0:
1612 for (i = 0; i < P - 2; ++i)
1613 {
1614 maparray[i] = GetMode(i + 2, 0, 0);
1615 }
1616 if (edgeOrient == eBackwards)
1617 {
1618 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1619 {
1620 signarray[i] = -1;
1621 }
1622 }
1623 break;
1624 case 1:
1625 for (i = 0; i < Q - 2; ++i)
1626 {
1627 maparray[i] = GetMode(1, i + 1, 0);
1628 }
1629 if (edgeOrient == eBackwards)
1630 {
1631 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1632 {
1633 signarray[i] = -1;
1634 }
1635 }
1636 break;
1637 case 2:
1638 for (i = 0; i < Q - 2; ++i)
1639 {
1640 maparray[i] = GetMode(0, i + 2, 0);
1641 }
1642 if (edgeOrient == eBackwards)
1643 {
1644 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1645 {
1646 signarray[i] = -1;
1647 }
1648 }
1649 break;
1650 case 3:
1651 for (i = 0; i < R - 2; ++i)
1652 {
1653 maparray[i] = GetMode(0, 0, i + 2);
1654 }
1655 if (edgeOrient == eBackwards)
1656 {
1657 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1658 {
1659 signarray[i] = -1;
1660 }
1661 }
1662 break;
1663 case 4:
1664 for (i = 0; i < R - 2; ++i)
1665 {
1666 maparray[i] = GetMode(1, 0, i + 1);
1667 }
1668 if (edgeOrient == eBackwards)
1669 {
1670 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1671 {
1672 signarray[i] = -1;
1673 }
1674 }
1675 break;
1676 case 5:
1677 for (i = 0; i < R - 2; ++i)
1678 {
1679 maparray[i] = GetMode(0, 1, i + 1);
1680 }
1681 if (edgeOrient == eBackwards)
1682 {
1683 for (i = 1; i < nEdgeIntCoeffs; i += 2)
1684 {
1685 signarray[i] = -1;
1686 }
1687 }
1688 break;
1689 default:
1690 ASSERTL0(false, "Edge not defined.");
1691 break;
1692 }
1693}
int v_GetEdgeNcoeffs(const int i) const override
Definition: StdTetExp.cpp:1080

References ASSERTL0, Nektar::StdRegions::eBackwards, GetMode(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LibUtilities::P, and v_GetEdgeNcoeffs().

◆ v_GetEdgeNcoeffs()

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

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 1080 of file StdTetExp.cpp.

1081{
1082 ASSERTL2((i >= 0) && (i <= 5), "edge id is out of range");
1083 int P = m_base[0]->GetNumModes();
1084 int Q = m_base[1]->GetNumModes();
1085 int R = m_base[2]->GetNumModes();
1086
1087 if (i == 0)
1088 {
1089 return P;
1090 }
1091 else if (i == 1 || i == 2)
1092 {
1093 return Q;
1094 }
1095 else
1096 {
1097 return R;
1098 }
1099}
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:265

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

Referenced by v_GetEdgeInteriorToElementMap().

◆ v_GetElmtTraceToTraceMap()

void Nektar::StdRegions::StdTetExp::v_GetElmtTraceToTraceMap ( const unsigned int  tid,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray,
Orientation  traceOrient = eForwards,
int  P = -1,
int  Q = -1 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1465 of file StdTetExp.cpp.

1469{
1470 int nummodesA = 0, nummodesB = 0, i, j, k, idx;
1471
1473 "Method only implemented for Modified_A BasisType (x "
1474 "direction), Modified_B BasisType (y direction), and "
1475 "Modified_C BasisType(z direction)");
1476
1477 int nFaceCoeffs = 0;
1478
1479 switch (fid)
1480 {
1481 case 0:
1482 nummodesA = m_base[0]->GetNumModes();
1483 nummodesB = m_base[1]->GetNumModes();
1484 break;
1485 case 1:
1486 nummodesA = m_base[0]->GetNumModes();
1487 nummodesB = m_base[2]->GetNumModes();
1488 break;
1489 case 2:
1490 case 3:
1491 nummodesA = m_base[1]->GetNumModes();
1492 nummodesB = m_base[2]->GetNumModes();
1493 break;
1494 default:
1495 ASSERTL0(false, "fid must be between 0 and 3");
1496 }
1497
1498 if (P == -1)
1499 {
1500 P = nummodesA;
1501 Q = nummodesB;
1502 }
1503
1504 nFaceCoeffs = P * (2 * Q - P + 1) / 2;
1505
1506 // Allocate the map array and sign array; set sign array to ones (+)
1507 if (maparray.size() != nFaceCoeffs)
1508 {
1509 maparray = Array<OneD, unsigned int>(nFaceCoeffs, 1);
1510 }
1511
1512 if (signarray.size() != nFaceCoeffs)
1513 {
1514 signarray = Array<OneD, int>(nFaceCoeffs, 1);
1515 }
1516 else
1517 {
1518 fill(signarray.get(), signarray.get() + nFaceCoeffs, 1);
1519 }
1520
1521 // zero signmap and set maparray to zero if elemental
1522 // modes are not as large as face modesl
1523 idx = 0;
1524 int cnt = 0;
1525 int minPA = min(nummodesA, P);
1526 int minQB = min(nummodesB, Q);
1527
1528 for (j = 0; j < minPA; ++j)
1529 {
1530 // set maparray
1531 for (k = 0; k < minQB - j; ++k, ++cnt)
1532 {
1533 maparray[idx++] = cnt;
1534 }
1535
1536 cnt += nummodesB - minQB;
1537
1538 for (k = nummodesB - j; k < Q - j; ++k)
1539 {
1540 signarray[idx] = 0.0;
1541 maparray[idx++] = maparray[0];
1542 }
1543 }
1544
1545 for (j = nummodesA; j < P; ++j)
1546 {
1547 for (k = 0; k < Q - j; ++k)
1548 {
1549 signarray[idx] = 0.0;
1550 maparray[idx++] = maparray[0];
1551 }
1552 }
1553
1554 if (faceOrient == eDir1BwdDir1_Dir2FwdDir2)
1555 {
1556 idx = 0;
1557 for (i = 0; i < P; ++i)
1558 {
1559 for (j = 0; j < Q - i; ++j, idx++)
1560 {
1561 if (i > 1)
1562 {
1563 signarray[idx] = (i % 2 ? -1 : 1);
1564 }
1565 }
1566 }
1567
1568 swap(maparray[0], maparray[Q]);
1569
1570 for (i = 1; i < Q - 1; ++i)
1571 {
1572 swap(maparray[i + 1], maparray[Q + i]);
1573 }
1574 }
1575}
bool v_IsBoundaryInteriorExpansion() const override
Definition: StdTetExp.cpp:1188

References ASSERTL0, ASSERTL1, Nektar::StdRegions::eDir1BwdDir1_Dir2FwdDir2, Nektar::StdRegions::StdExpansion::m_base, Nektar::LibUtilities::P, and v_IsBoundaryInteriorExpansion().

◆ v_GetInteriorMap()

void Nektar::StdRegions::StdTetExp::v_GetInteriorMap ( Array< OneD, unsigned int > &  outarray)
overrideprotectedvirtual

Maps interior modes of an edge to the elemental modes. List of all interior modes in the expansion.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1279 of file StdTetExp.cpp.

1280{
1283 "BasisType is not a boundary interior form");
1286 "BasisType is not a boundary interior form");
1289 "BasisType is not a boundary interior form");
1290
1291 int P = m_base[0]->GetNumModes();
1292 int Q = m_base[1]->GetNumModes();
1293 int R = m_base[2]->GetNumModes();
1294
1295 int nIntCoeffs = m_ncoeffs - NumBndryCoeffs();
1296
1297 if (outarray.size() != nIntCoeffs)
1298 {
1299 outarray = Array<OneD, unsigned int>(nIntCoeffs);
1300 }
1301
1302 int idx = 0;
1303 for (int i = 2; i < P; ++i)
1304 {
1305 for (int j = 1; j < Q - i - 1; ++j)
1306 {
1307 for (int k = 1; k < R - i - j; ++k)
1308 {
1309 outarray[idx++] = GetMode(i, j, k);
1310 }
1311 }
1312 }
1313}

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

◆ v_GetNedges()

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

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 960 of file StdTetExp.cpp.

961{
962 return 6;
963}

◆ v_GetNtraces()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 965 of file StdTetExp.cpp.

966{
967 return 4;
968}

◆ v_GetNverts()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 955 of file StdTetExp.cpp.

956{
957 return 4;
958}

◆ v_GetSimplexEquiSpacedConnectivity()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2233 of file StdTetExp.cpp.

2235{
2236 int np0 = m_base[0]->GetNumPoints();
2237 int np1 = m_base[1]->GetNumPoints();
2238 int np2 = m_base[2]->GetNumPoints();
2239 int np = max(np0, max(np1, np2));
2240
2241 conn = Array<OneD, int>(4 * (np - 1) * (np - 1) * (np - 1));
2242
2243 int row = 0;
2244 int rowp1 = 0;
2245 int plane = 0;
2246 int row1 = 0;
2247 int row1p1 = 0;
2248 int planep1 = 0;
2249 int cnt = 0;
2250 for (int i = 0; i < np - 1; ++i)
2251 {
2252 planep1 += (np - i) * (np - i + 1) / 2;
2253 row = 0; // current plane row offset
2254 rowp1 = 0; // current plane row plus one offset
2255 row1 = 0; // next plane row offset
2256 row1p1 = 0; // nex plane row plus one offset
2257 for (int j = 0; j < np - i - 1; ++j)
2258 {
2259 rowp1 += np - i - j;
2260 row1p1 += np - i - j - 1;
2261 for (int k = 0; k < np - i - j - 2; ++k)
2262 {
2263 conn[cnt++] = plane + row + k + 1;
2264 conn[cnt++] = plane + row + k;
2265 conn[cnt++] = plane + rowp1 + k;
2266 conn[cnt++] = planep1 + row1 + k;
2267
2268 conn[cnt++] = plane + row + k + 1;
2269 conn[cnt++] = plane + rowp1 + k + 1;
2270 conn[cnt++] = planep1 + row1 + k + 1;
2271 conn[cnt++] = planep1 + row1 + k;
2272
2273 conn[cnt++] = plane + rowp1 + k + 1;
2274 conn[cnt++] = plane + row + k + 1;
2275 conn[cnt++] = plane + rowp1 + k;
2276 conn[cnt++] = planep1 + row1 + k;
2277
2278 conn[cnt++] = planep1 + row1 + k;
2279 conn[cnt++] = planep1 + row1p1 + k;
2280 conn[cnt++] = plane + rowp1 + k;
2281 conn[cnt++] = plane + rowp1 + k + 1;
2282
2283 conn[cnt++] = planep1 + row1 + k;
2284 conn[cnt++] = planep1 + row1p1 + k;
2285 conn[cnt++] = planep1 + row1 + k + 1;
2286 conn[cnt++] = plane + rowp1 + k + 1;
2287
2288 if (k < np - i - j - 3)
2289 {
2290 conn[cnt++] = plane + rowp1 + k + 1;
2291 conn[cnt++] = planep1 + row1p1 + k + 1;
2292 conn[cnt++] = planep1 + row1 + k + 1;
2293 conn[cnt++] = planep1 + row1p1 + k;
2294 }
2295 }
2296
2297 conn[cnt++] = plane + row + np - i - j - 1;
2298 conn[cnt++] = plane + row + np - i - j - 2;
2299 conn[cnt++] = plane + rowp1 + np - i - j - 2;
2300 conn[cnt++] = planep1 + row1 + np - i - j - 2;
2301
2302 row += np - i - j;
2303 row1 += np - i - j - 1;
2304 }
2305 plane += (np - i) * (np - i + 1) / 2;
2306 }
2307}

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

◆ v_GetTraceBasisKey()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1132 of file StdTetExp.cpp.

1134{
1135 ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
1136 ASSERTL2(k == 0 || k == 1, "face direction out of range");
1137
1138 int dir = k;
1139 switch (i)
1140 {
1141 case 0:
1142 dir = k;
1143 break;
1144 case 1:
1145 dir = 2 * k;
1146 break;
1147 case 2:
1148 case 3:
1149 dir = k + 1;
1150 break;
1151 }
1152
1154 m_base[dir]->GetNumPoints(),
1155 m_base[dir]->GetNumModes());
1156}
LibUtilities::BasisKey EvaluateTriFaceBasisKey(const int facedir, const LibUtilities::BasisType faceDirBasisType, const int numpoints, const int nummodes)

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

◆ v_GetTraceCoeffMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1373 of file StdTetExp.cpp.

1375{
1376 int i, j, k;
1377 int P = 0, Q = 0, idx = 0;
1378 int nFaceCoeffs = 0;
1379
1380 switch (fid)
1381 {
1382 case 0:
1383 P = m_base[0]->GetNumModes();
1384 Q = m_base[1]->GetNumModes();
1385 break;
1386 case 1:
1387 P = m_base[0]->GetNumModes();
1388 Q = m_base[2]->GetNumModes();
1389 break;
1390 case 2:
1391 case 3:
1392 P = m_base[1]->GetNumModes();
1393 Q = m_base[2]->GetNumModes();
1394 break;
1395 default:
1396 ASSERTL0(false, "fid must be between 0 and 3");
1397 }
1398
1399 nFaceCoeffs = P * (2 * Q - P + 1) / 2;
1400
1401 if (maparray.size() != nFaceCoeffs)
1402 {
1403 maparray = Array<OneD, unsigned int>(nFaceCoeffs);
1404 }
1405
1406 switch (fid)
1407 {
1408 case 0:
1409 idx = 0;
1410 for (i = 0; i < P; ++i)
1411 {
1412 for (j = 0; j < Q - i; ++j)
1413 {
1414 maparray[idx++] = GetMode(i, j, 0);
1415 }
1416 }
1417 break;
1418 case 1:
1419 idx = 0;
1420 for (i = 0; i < P; ++i)
1421 {
1422 for (k = 0; k < Q - i; ++k)
1423 {
1424 maparray[idx++] = GetMode(i, 0, k);
1425 }
1426 }
1427 break;
1428 case 2:
1429 idx = 0;
1430 for (j = 0; j < P - 1; ++j)
1431 {
1432 for (k = 0; k < Q - 1 - j; ++k)
1433 {
1434 maparray[idx++] = GetMode(1, j, k);
1435 // Incorporate modes from zeroth plane where needed.
1436 if (j == 0 && k == 0)
1437 {
1438 maparray[idx++] = GetMode(0, 0, 1);
1439 }
1440 if (j == 0 && k == Q - 2)
1441 {
1442 for (int r = 0; r < Q - 1; ++r)
1443 {
1444 maparray[idx++] = GetMode(0, 1, r);
1445 }
1446 }
1447 }
1448 }
1449 break;
1450 case 3:
1451 idx = 0;
1452 for (j = 0; j < P; ++j)
1453 {
1454 for (k = 0; k < Q - j; ++k)
1455 {
1456 maparray[idx++] = GetMode(0, j, k);
1457 }
1458 }
1459 break;
1460 default:
1461 ASSERTL0(false, "Element map not available.");
1462 }
1463}

References ASSERTL0, GetMode(), Nektar::StdRegions::StdExpansion::m_base, and Nektar::LibUtilities::P.

◆ v_GetTraceInteriorToElementMap()

void Nektar::StdRegions::StdTetExp::v_GetTraceInteriorToElementMap ( const int  tid,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray,
const Orientation  traceOrient = eDir1FwdDir1_Dir2FwdDir2 
)
overrideprotectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1695 of file StdTetExp.cpp.

1698{
1699 int i, j, idx, k;
1700 const int P = m_base[0]->GetNumModes();
1701 const int Q = m_base[1]->GetNumModes();
1702 const int R = m_base[2]->GetNumModes();
1703
1704 const int nFaceIntCoeffs = v_GetTraceIntNcoeffs(fid);
1705
1706 if (maparray.size() != nFaceIntCoeffs)
1707 {
1708 maparray = Array<OneD, unsigned int>(nFaceIntCoeffs);
1709 }
1710
1711 if (signarray.size() != nFaceIntCoeffs)
1712 {
1713 signarray = Array<OneD, int>(nFaceIntCoeffs, 1);
1714 }
1715 else
1716 {
1717 fill(signarray.get(), signarray.get() + nFaceIntCoeffs, 1);
1718 }
1719
1720 switch (fid)
1721 {
1722 case 0:
1723 idx = 0;
1724 for (i = 2; i < P; ++i)
1725 {
1726 for (j = 1; j < Q - i; ++j)
1727 {
1728 if ((int)faceOrient == 7)
1729 {
1730 signarray[idx] = (i % 2 ? -1 : 1);
1731 }
1732 maparray[idx++] = GetMode(i, j, 0);
1733 }
1734 }
1735 break;
1736 case 1:
1737 idx = 0;
1738 for (i = 2; i < P; ++i)
1739 {
1740 for (k = 1; k < R - i; ++k)
1741 {
1742 if ((int)faceOrient == 7)
1743 {
1744 signarray[idx] = (i % 2 ? -1 : 1);
1745 }
1746 maparray[idx++] = GetMode(i, 0, k);
1747 }
1748 }
1749 break;
1750 case 2:
1751 idx = 0;
1752 for (j = 1; j < Q - 1; ++j)
1753 {
1754 for (k = 1; k < R - 1 - j; ++k)
1755 {
1756 if ((int)faceOrient == 7)
1757 {
1758 signarray[idx] = ((j + 1) % 2 ? -1 : 1);
1759 }
1760 maparray[idx++] = GetMode(1, j, k);
1761 }
1762 }
1763 break;
1764 case 3:
1765 idx = 0;
1766 for (j = 2; j < Q; ++j)
1767 {
1768 for (k = 1; k < R - j; ++k)
1769 {
1770 if ((int)faceOrient == 7)
1771 {
1772 signarray[idx] = (j % 2 ? -1 : 1);
1773 }
1774 maparray[idx++] = GetMode(0, j, k);
1775 }
1776 }
1777 break;
1778 default:
1779 ASSERTL0(false, "Face interior map not available.");
1780 break;
1781 }
1782}
int v_GetTraceIntNcoeffs(const int i) const override
Definition: StdTetExp.cpp:1041

References ASSERTL0, GetMode(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LibUtilities::P, and v_GetTraceIntNcoeffs().

◆ v_GetTraceIntNcoeffs()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 1041 of file StdTetExp.cpp.

1042{
1043 ASSERTL2((i >= 0) && (i <= 3), "face id is out of range");
1044 int Pi = m_base[0]->GetNumModes() - 2;
1045 int Qi = m_base[1]->GetNumModes() - 2;
1046 int Ri = m_base[2]->GetNumModes() - 2;
1047
1048 if ((i == 0))
1049 {
1050 return Pi * (2 * Qi - Pi - 1) / 2;
1051 }
1052 else if ((i == 1))
1053 {
1054 return Pi * (2 * Ri - Pi - 1) / 2;
1055 }
1056 else
1057 {
1058 return Qi * (2 * Ri - Qi - 1) / 2;
1059 }
1060}

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

Referenced by v_GetTraceInteriorToElementMap().

◆ v_GetTraceNcoeffs()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 1015 of file StdTetExp.cpp.

1016{
1017 ASSERTL2((i >= 0) && (i <= 3), "face id is out of range");
1018 int nFaceCoeffs = 0;
1019 int nummodesA, nummodesB, P, Q;
1020 if (i == 0)
1021 {
1022 nummodesA = GetBasisNumModes(0);
1023 nummodesB = GetBasisNumModes(1);
1024 }
1025 else if ((i == 1) || (i == 2))
1026 {
1027 nummodesA = GetBasisNumModes(0);
1028 nummodesB = GetBasisNumModes(2);
1029 }
1030 else
1031 {
1032 nummodesA = GetBasisNumModes(1);
1033 nummodesB = GetBasisNumModes(2);
1034 }
1035 P = nummodesA - 1;
1036 Q = nummodesB - 1;
1037 nFaceCoeffs = Q + 1 + (P * (1 + 2 * Q - P)) / 2;
1038 return nFaceCoeffs;
1039}
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:169

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

◆ v_GetTraceNumModes()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 921 of file StdTetExp.cpp.

924{
925 int nummodes[3] = {m_base[0]->GetNumModes(), m_base[1]->GetNumModes(),
926 m_base[2]->GetNumModes()};
927 switch (fid)
928 {
929 case 0:
930 {
931 numModes0 = nummodes[0];
932 numModes1 = nummodes[1];
933 }
934 break;
935 case 1:
936 {
937 numModes0 = nummodes[0];
938 numModes1 = nummodes[2];
939 }
940 break;
941 case 2:
942 case 3:
943 {
944 numModes0 = nummodes[1];
945 numModes1 = nummodes[2];
946 }
947 break;
948 }
949}

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

◆ v_GetTraceNumPoints()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 1062 of file StdTetExp.cpp.

1063{
1064 ASSERTL2(i >= 0 && i <= 3, "face id is out of range");
1065
1066 if (i == 0)
1067 {
1068 return m_base[0]->GetNumPoints() * m_base[1]->GetNumPoints();
1069 }
1070 else if (i == 1)
1071 {
1072 return m_base[0]->GetNumPoints() * m_base[2]->GetNumPoints();
1073 }
1074 else
1075 {
1076 return m_base[1]->GetNumPoints() * m_base[2]->GetNumPoints();
1077 }
1078}

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

◆ v_GetTracePointsKey()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1101 of file StdTetExp.cpp.

1103{
1104 ASSERTL2(i >= 0 && i <= 3, "face id is out of range");
1105 ASSERTL2(j == 0 || j == 1, "face direction is out of range");
1106
1107 if (i == 0)
1108 {
1109 return m_base[j]->GetPointsKey();
1110 }
1111 else if (i == 1)
1112 {
1113 return m_base[2 * j]->GetPointsKey();
1114 }
1115 else
1116 {
1117 return m_base[j + 1]->GetPointsKey();
1118 }
1119}

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

◆ v_GetVertexMap()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1198 of file StdTetExp.cpp.

1199{
1203 "Mapping not defined for this type of basis");
1204
1205 int localDOF = 0;
1206 if (useCoeffPacking == true) // follow packing of coefficients i.e q,r,p
1207 {
1208 switch (localVertexId)
1209 {
1210 case 0:
1211 {
1212 localDOF = GetMode(0, 0, 0);
1213 break;
1214 }
1215 case 1:
1216 {
1217 localDOF = GetMode(0, 0, 1);
1218 break;
1219 }
1220 case 2:
1221 {
1222 localDOF = GetMode(0, 1, 0);
1223 break;
1224 }
1225 case 3:
1226 {
1227 localDOF = GetMode(1, 0, 0);
1228 break;
1229 }
1230 default:
1231 {
1232 ASSERTL0(false, "Vertex ID must be between 0 and 3");
1233 break;
1234 }
1235 }
1236 }
1237 else
1238 {
1239 switch (localVertexId)
1240 {
1241 case 0:
1242 {
1243 localDOF = GetMode(0, 0, 0);
1244 break;
1245 }
1246 case 1:
1247 {
1248 localDOF = GetMode(1, 0, 0);
1249 break;
1250 }
1251 case 2:
1252 {
1253 localDOF = GetMode(0, 1, 0);
1254 break;
1255 }
1256 case 3:
1257 {
1258 localDOF = GetMode(0, 0, 1);
1259 break;
1260 }
1261 default:
1262 {
1263 ASSERTL0(false, "Vertex ID must be between 0 and 3");
1264 break;
1265 }
1266 }
1267 }
1268
1269 return localDOF;
1270}

References ASSERTL0, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::LibUtilities::eModified_C, Nektar::StdRegions::StdExpansion::GetBasisType(), and GetMode().

◆ v_IProductWRTBase()

void Nektar::StdRegions::StdTetExp::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} (\eta_{1i}) \psi_{pq}^{b} (\eta_{2j}) \psi_{pqr}^{c} (\eta_{3k}) w_i w_j w_k u(\eta_{1,i} \eta_{2,j} \eta_{3,k}) J_{i,j,k}\\ & = & \sum_{i=0}^{nq_0} \psi_p^a(\eta_{1,i}) \sum_{j=0}^{nq_1} \psi_{pq}^b(\eta_{2,j}) \sum_{k=0}^{nq_2} \psi_{pqr}^c u(\eta_{1i},\eta_{2j},\eta_{3k}) J_{i,j,k} \end{array} \)
where

\( \phi_{pqr} (\xi_1 , \xi_2 , \xi_3) = \psi_p^a (\eta_1) \psi_{pq}^b (\eta_2) \psi_{pqr}^c (\eta_3) \)

which can be implemented as
\(f_{pqr} (\xi_{3k}) = \sum_{k=0}^{nq_3} \psi_{pqr}^c u(\eta_{1i},\eta_{2j},\eta_{3k}) J_{i,j,k} = {\bf B_3 U} \)
\( g_{pq} (\xi_{3k}) = \sum_{j=0}^{nq_1} \psi_{pq}^b (\xi_{2j}) f_{pqr} (\xi_{3k}) = {\bf B_2 F} \)
\( (\phi_{pqr}, u)_{\delta} = \sum_{k=0}^{nq_0} \psi_{p}^a (\xi_{3k}) g_{pq} (\xi_{3k}) = {\bf B_1 G} \)

Parameters
inarrayFunction evaluated at physical collocation points.
outarrayInner product with respect to each basis function over the element.

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TetExp.

Definition at line 466 of file StdTetExp.cpp.

468{
471 "Basis[1] is not a general tensor type");
472
475 "Basis[2] is not a general tensor type");
476
477 if (m_base[0]->Collocation() && m_base[1]->Collocation())
478 {
479 MultiplyByQuadratureMetric(inarray, outarray);
480 }
481 else
482 {
483 StdTetExp::v_IProductWRTBase_SumFac(inarray, outarray);
484 }
485}
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:723
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
Definition: StdTetExp.cpp:493

References ASSERTL1, Nektar::LibUtilities::eModified_B, Nektar::LibUtilities::eModified_C, Nektar::LibUtilities::eOrtho_B, Nektar::LibUtilities::eOrtho_C, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::MultiplyByQuadratureMetric(), and v_IProductWRTBase_SumFac().

Referenced by v_FwdTrans().

◆ v_IProductWRTBase_SumFac()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TetExp.

Definition at line 493 of file StdTetExp.cpp.

496{
497 int nquad0 = m_base[0]->GetNumPoints();
498 int nquad1 = m_base[1]->GetNumPoints();
499 int nquad2 = m_base[2]->GetNumPoints();
500 int order0 = m_base[0]->GetNumModes();
501 int order1 = m_base[1]->GetNumModes();
502
503 Array<OneD, NekDouble> wsp(nquad1 * nquad2 * order0 +
504 nquad2 * order0 * (2 * order1 - order0 + 1) / 2);
505
506 if (multiplybyweights)
507 {
508 Array<OneD, NekDouble> tmp(nquad0 * nquad1 * nquad2);
509 MultiplyByQuadratureMetric(inarray, tmp);
510
512 m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
513 tmp, outarray, wsp, true, true, true);
514 }
515 else
516 {
518 m_base[0]->GetBdata(), m_base[1]->GetBdata(), m_base[2]->GetBdata(),
519 inarray, outarray, wsp, true, true, true);
520 }
521}
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(), and Nektar::StdRegions::StdNodalTetExp::v_IProductWRTBase_SumFac().

◆ v_IProductWRTBase_SumFacKernel()

void Nektar::StdRegions::StdTetExp::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

Implements Nektar::StdRegions::StdExpansion3D.

Definition at line 523 of file StdTetExp.cpp.

532{
533 int nquad0 = m_base[0]->GetNumPoints();
534 int nquad1 = m_base[1]->GetNumPoints();
535 int nquad2 = m_base[2]->GetNumPoints();
536
537 int order0 = m_base[0]->GetNumModes();
538 int order1 = m_base[1]->GetNumModes();
539 int order2 = m_base[2]->GetNumModes();
540
541 Array<OneD, NekDouble> tmp1 = wsp;
542 Array<OneD, NekDouble> tmp2 = wsp + nquad1 * nquad2 * order0;
543
544 int i, j, mode, mode1, cnt;
545
546 // Inner product with respect to the '0' direction
547 Blas::Dgemm('T', 'N', nquad1 * nquad2, order0, nquad0, 1.0, inarray.get(),
548 nquad0, base0.get(), nquad0, 0.0, tmp1.get(), nquad1 * nquad2);
549
550 // Inner product with respect to the '1' direction
551 for (mode = i = 0; i < order0; ++i)
552 {
553 Blas::Dgemm('T', 'N', nquad2, order1 - i, nquad1, 1.0,
554 tmp1.get() + i * nquad1 * nquad2, nquad1,
555 base1.get() + mode * nquad1, nquad1, 0.0,
556 tmp2.get() + mode * nquad2, nquad2);
557 mode += order1 - i;
558 }
559
560 // fix for modified basis for base singular vertex
562 {
563 // base singular vertex and singular edge (1+b)/2
564 //(1+a)/2 components (makes tmp[nquad2] entry into (1+b)/2)
565 Blas::Dgemv('T', nquad1, nquad2, 1.0, tmp1.get() + nquad1 * nquad2,
566 nquad1, base1.get() + nquad1, 1, 1.0, tmp2.get() + nquad2,
567 1);
568 }
569
570 // Inner product with respect to the '2' direction
571 mode = mode1 = cnt = 0;
572 for (i = 0; i < order0; ++i)
573 {
574 for (j = 0; j < order1 - i; ++j, ++cnt)
575 {
576 Blas::Dgemv('T', nquad2, order2 - i - j, 1.0,
577 base2.get() + mode * nquad2, nquad2,
578 tmp2.get() + cnt * nquad2, 1, 0.0,
579 outarray.get() + mode1, 1);
580 mode += order2 - i - j;
581 mode1 += order2 - i - j;
582 }
583 // increment mode in case order1!=order2
584 for (j = order1 - i; j < order2 - i; ++j)
585 {
586 mode += order2 - i - j;
587 }
588 }
589
590 // fix for modified basis for top singular vertex component
591 // Already have evaluated (1+c)/2 (1-b)/2 (1-a)/2
593 {
594 // add in (1+c)/2 (1+b)/2 component
595 outarray[1] +=
596 Blas::Ddot(nquad2, base2.get() + nquad2, 1, &tmp2[nquad2], 1);
597
598 // add in (1+c)/2 (1-b)/2 (1+a)/2 component
599 outarray[1] += Blas::Ddot(nquad2, base2.get() + nquad2, 1,
600 &tmp2[nquad2 * order1], 1);
601 }
602}
static double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = .
Definition: Blas.hpp:163

References Blas::Ddot(), Blas::Dgemm(), Blas::Dgemv(), Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), and Nektar::StdRegions::StdExpansion::m_base.

◆ v_IProductWRTDerivBase()

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

Definition at line 604 of file StdTetExp.cpp.

607{
608 StdTetExp::v_IProductWRTDerivBase_SumFac(dir, inarray, outarray);
609}
void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdTetExp.cpp:617

References v_IProductWRTDerivBase_SumFac().

◆ v_IProductWRTDerivBase_SumFac()

void Nektar::StdRegions::StdTetExp::v_IProductWRTDerivBase_SumFac ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual
Parameters
inarrayFunction evaluated at physical collocation points.
outarrayInner product with respect to each basis function over the element.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 617 of file StdTetExp.cpp.

620{
621 int i;
622 int nquad0 = m_base[0]->GetNumPoints();
623 int nquad1 = m_base[1]->GetNumPoints();
624 int nquad2 = m_base[2]->GetNumPoints();
625 int nqtot = nquad0 * nquad1 * nquad2;
626 int nmodes0 = m_base[0]->GetNumModes();
627 int nmodes1 = m_base[1]->GetNumModes();
628 int wspsize = nquad0 + nquad1 + nquad2 + max(nqtot, m_ncoeffs) +
629 nquad1 * nquad2 * nmodes0 +
630 nquad2 * nmodes0 * (2 * nmodes1 - nmodes0 + 1) / 2;
631
632 Array<OneD, NekDouble> gfac0(wspsize);
633 Array<OneD, NekDouble> gfac1(gfac0 + nquad0);
634 Array<OneD, NekDouble> gfac2(gfac1 + nquad1);
635 Array<OneD, NekDouble> tmp0(gfac2 + nquad2);
636 Array<OneD, NekDouble> wsp(tmp0 + max(nqtot, m_ncoeffs));
637
638 const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
639 const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
640 const Array<OneD, const NekDouble> &z2 = m_base[2]->GetZ();
641
642 // set up geometric factor: (1+z0)/2
643 for (i = 0; i < nquad0; ++i)
644 {
645 gfac0[i] = 0.5 * (1 + z0[i]);
646 }
647
648 // set up geometric factor: 2/(1-z1)
649 for (i = 0; i < nquad1; ++i)
650 {
651 gfac1[i] = 2.0 / (1 - z1[i]);
652 }
653
654 // Set up geometric factor: 2/(1-z2)
655 for (i = 0; i < nquad2; ++i)
656 {
657 gfac2[i] = 2.0 / (1 - z2[i]);
658 }
659
660 // Derivative in first direction is always scaled as follows
661 for (i = 0; i < nquad1 * nquad2; ++i)
662 {
663 Vmath::Smul(nquad0, gfac1[i % nquad1], &inarray[0] + i * nquad0, 1,
664 &tmp0[0] + i * nquad0, 1);
665 }
666 for (i = 0; i < nquad2; ++i)
667 {
668 Vmath::Smul(nquad0 * nquad1, gfac2[i], &tmp0[0] + i * nquad0 * nquad1,
669 1, &tmp0[0] + i * nquad0 * nquad1, 1);
670 }
671
672 MultiplyByQuadratureMetric(tmp0, tmp0);
673
674 switch (dir)
675 {
676 case 0:
677 {
679 m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
680 m_base[2]->GetBdata(), tmp0, outarray, wsp, false, true, true);
681 }
682 break;
683 case 1:
684 {
685 Array<OneD, NekDouble> tmp3(m_ncoeffs);
686
687 for (i = 0; i < nquad1 * nquad2; ++i)
688 {
689 Vmath::Vmul(nquad0, &gfac0[0], 1, &tmp0[0] + i * nquad0, 1,
690 &tmp0[0] + i * nquad0, 1);
691 }
692
694 m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
695 m_base[2]->GetBdata(), tmp0, tmp3, wsp, false, true, true);
696
697 for (i = 0; i < nquad2; ++i)
698 {
699 Vmath::Smul(nquad0 * nquad1, gfac2[i],
700 &inarray[0] + i * nquad0 * nquad1, 1,
701 &tmp0[0] + i * nquad0 * nquad1, 1);
702 }
703 MultiplyByQuadratureMetric(tmp0, tmp0);
705 m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
706 m_base[2]->GetBdata(), tmp0, outarray, wsp, true, false, true);
707 Vmath::Vadd(m_ncoeffs, &tmp3[0], 1, &outarray[0], 1, &outarray[0],
708 1);
709 }
710 break;
711 case 2:
712 {
713 Array<OneD, NekDouble> tmp3(m_ncoeffs);
714 Array<OneD, NekDouble> tmp4(m_ncoeffs);
715 for (i = 0; i < nquad1; ++i)
716 {
717 gfac1[i] = (1 + z1[i]) / 2;
718 }
719
720 for (i = 0; i < nquad1 * nquad2; ++i)
721 {
722 Vmath::Vmul(nquad0, &gfac0[0], 1, &tmp0[0] + i * nquad0, 1,
723 &tmp0[0] + i * nquad0, 1);
724 }
726 m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
727 m_base[2]->GetBdata(), tmp0, tmp3, wsp, false, true, true);
728
729 for (i = 0; i < nquad2; ++i)
730 {
731 Vmath::Smul(nquad0 * nquad1, gfac2[i],
732 &inarray[0] + i * nquad0 * nquad1, 1,
733 &tmp0[0] + i * nquad0 * nquad1, 1);
734 }
735 for (i = 0; i < nquad1 * nquad2; ++i)
736 {
737 Vmath::Smul(nquad0, gfac1[i % nquad1], &tmp0[0] + i * nquad0, 1,
738 &tmp0[0] + i * nquad0, 1);
739 }
740 MultiplyByQuadratureMetric(tmp0, tmp0);
742 m_base[0]->GetBdata(), m_base[1]->GetDbdata(),
743 m_base[2]->GetBdata(), tmp0, tmp4, wsp, true, false, true);
744
745 MultiplyByQuadratureMetric(inarray, tmp0);
747 m_base[0]->GetBdata(), m_base[1]->GetBdata(),
748 m_base[2]->GetDbdata(), tmp0, outarray, wsp, true, true, false);
749
750 Vmath::Vadd(m_ncoeffs, &tmp3[0], 1, &outarray[0], 1, &outarray[0],
751 1);
752 Vmath::Vadd(m_ncoeffs, &tmp4[0], 1, &outarray[0], 1, &outarray[0],
753 1);
754 }
755 break;
756 default:
757 {
758 ASSERTL1(false, "input dir is out of range");
759 }
760 break;
761 }
762}
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.hpp:72
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180

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

Referenced by v_IProductWRTDerivBase(), and Nektar::StdRegions::StdNodalTetExp::v_IProductWRTDerivBase_SumFac().

◆ v_IsBoundaryInteriorExpansion()

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

◆ v_LocCollapsedToLocCoord()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 800 of file StdTetExp.cpp.

802{
803 xi[1] = (1.0 + eta[0]) * (1.0 - eta[2]) * 0.5 - 1.0;
804 xi[0] = (1.0 + eta[0]) * (-xi[1] - eta[2]) * 0.5 - 1.0;
805 xi[2] = eta[2];
806}

◆ v_LocCoordToLocCollapsed()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 768 of file StdTetExp.cpp.

770{
771 NekDouble d2 = 1.0 - xi[2];
772 NekDouble d12 = -xi[1] - xi[2];
773 if (fabs(d2) < NekConstants::kNekZeroTol)
774 {
775 if (d2 >= 0.)
776 {
778 }
779 else
780 {
782 }
783 }
784 if (fabs(d12) < NekConstants::kNekZeroTol)
785 {
786 if (d12 >= 0.)
787 {
789 }
790 else
791 {
793 }
794 }
795 eta[0] = 2.0 * (1.0 + xi[0]) / d12 - 1.0;
796 eta[1] = 2.0 * (1.0 + xi[1]) / d2 - 1.0;
797 eta[2] = xi[2];
798}
static const NekDouble kNekZeroTol

References Nektar::NekConstants::kNekZeroTol.

◆ v_MultiplyByStdQuadratureMetric()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1950 of file StdTetExp.cpp.

1953{
1954 int i, j;
1955
1956 int nquad0 = m_base[0]->GetNumPoints();
1957 int nquad1 = m_base[1]->GetNumPoints();
1958 int nquad2 = m_base[2]->GetNumPoints();
1959
1960 const Array<OneD, const NekDouble> &w0 = m_base[0]->GetW();
1961 const Array<OneD, const NekDouble> &w1 = m_base[1]->GetW();
1962 const Array<OneD, const NekDouble> &w2 = m_base[2]->GetW();
1963
1964 const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
1965 const Array<OneD, const NekDouble> &z2 = m_base[2]->GetZ();
1966
1967 // multiply by integration constants
1968 for (i = 0; i < nquad1 * nquad2; ++i)
1969 {
1970 Vmath::Vmul(nquad0, (NekDouble *)&inarray[0] + i * nquad0, 1, w0.get(),
1971 1, &outarray[0] + i * nquad0, 1);
1972 }
1973
1974 switch (m_base[1]->GetPointsType())
1975 {
1976 // (1,0) Jacobi Inner product.
1977 case LibUtilities::eGaussRadauMAlpha1Beta0:
1978 for (j = 0; j < nquad2; ++j)
1979 {
1980 for (i = 0; i < nquad1; ++i)
1981 {
1982 Blas::Dscal(nquad0, 0.5 * w1[i],
1983 &outarray[0] + i * nquad0 + j * nquad0 * nquad1,
1984 1);
1985 }
1986 }
1987 break;
1988
1989 default:
1990 for (j = 0; j < nquad2; ++j)
1991 {
1992 for (i = 0; i < nquad1; ++i)
1993 {
1994 Blas::Dscal(nquad0, 0.5 * (1 - z1[i]) * w1[i],
1995 &outarray[0] + i * nquad0 + j * nquad0 * nquad1,
1996 1);
1997 }
1998 }
1999 break;
2000 }
2001
2002 switch (m_base[2]->GetPointsType())
2003 {
2004 // (2,0) Jacobi inner product.
2005 case LibUtilities::eGaussRadauMAlpha2Beta0:
2006 for (i = 0; i < nquad2; ++i)
2007 {
2008 Blas::Dscal(nquad0 * nquad1, 0.25 * w2[i],
2009 &outarray[0] + i * nquad0 * nquad1, 1);
2010 }
2011 break;
2012 // (1,0) Jacobi inner product.
2013 case LibUtilities::eGaussRadauMAlpha1Beta0:
2014 for (i = 0; i < nquad2; ++i)
2015 {
2016 Blas::Dscal(nquad0 * nquad1, 0.25 * (1 - z2[i]) * w2[i],
2017 &outarray[0] + i * nquad0 * nquad1, 1);
2018 }
2019 break;
2020 default:
2021 for (i = 0; i < nquad2; ++i)
2022 {
2023 Blas::Dscal(nquad0 * nquad1,
2024 0.25 * (1 - z2[i]) * (1 - z2[i]) * w2[i],
2025 &outarray[0] + i * nquad0 * nquad1, 1);
2026 }
2027 break;
2028 }
2029}
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:205
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition: Blas.hpp:149

References Blas::Dscal(), Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::StdRegions::StdExpansion::m_base, and Vmath::Vmul().

◆ v_NumBndryCoeffs()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 975 of file StdTetExp.cpp.

976{
979 "BasisType is not a boundary interior form");
982 "BasisType is not a boundary interior form");
985 "BasisType is not a boundary interior form");
986
987 int P = m_base[0]->GetNumModes();
988 int Q = m_base[1]->GetNumModes();
989 int R = m_base[2]->GetNumModes();
990
992}
int getNumberOfBndCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:210

References ASSERTL1, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::LibUtilities::eModified_C, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::LibUtilities::StdTetData::getNumberOfBndCoefficients(), Nektar::StdRegions::StdExpansion::m_base, and Nektar::LibUtilities::P.

◆ v_NumDGBndryCoeffs()

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 994 of file StdTetExp.cpp.

995{
998 "BasisType is not a boundary interior form");
1001 "BasisType is not a boundary interior form");
1004 "BasisType is not a boundary interior form");
1005
1006 int P = m_base[0]->GetNumModes() - 1;
1007 int Q = m_base[1]->GetNumModes() - 1;
1008 int R = m_base[2]->GetNumModes() - 1;
1009
1010 return (Q + 1) + P * (1 + 2 * Q - P) / 2 // base face
1011 + (R + 1) + P * (1 + 2 * R - P) / 2 // front face
1012 + 2 * (R + 1) + Q * (1 + 2 * R - Q); // back two faces
1013}

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

◆ v_PhysDeriv() [1/2]

void Nektar::StdRegions::StdTetExp::v_PhysDeriv ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_dxi0,
Array< OneD, NekDouble > &  out_dxi1,
Array< OneD, NekDouble > &  out_dxi2 
)
overrideprotectedvirtual

Calculate the derivative of the physical points.

The derivative is evaluated at the nodal physical points. Derivatives with respect to the local Cartesian coordinates

\(\begin{Bmatrix} \frac {\partial} {\partial \xi_1} \\ \frac {\partial} {\partial \xi_2} \\ \frac {\partial} {\partial \xi_3} \end{Bmatrix} = \begin{Bmatrix} \frac 4 {(1-\eta_2)(1-\eta_3)} \frac \partial {\partial \eta_1} \ \ \frac {2(1+\eta_1)} {(1-\eta_2)(1-\eta_3)} \frac \partial {\partial \eta_1} + \frac 2 {1-\eta_3} \frac \partial {\partial \eta_3} \\ \frac {2(1 + \eta_1)} {2(1 - \eta_2)(1-\eta_3)} \frac \partial {\partial \eta_1} + \frac {1 + \eta_2} {1 - \eta_3} \frac \partial {\partial \eta_2} + \frac \partial {\partial \eta_3} \end{Bmatrix}\)

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TetExp.

Definition at line 83 of file StdTetExp.cpp.

87{
88 int Q0 = m_base[0]->GetNumPoints();
89 int Q1 = m_base[1]->GetNumPoints();
90 int Q2 = m_base[2]->GetNumPoints();
91 int Qtot = Q0 * Q1 * Q2;
92
93 // Compute the physical derivative
94 Array<OneD, NekDouble> out_dEta0(3 * Qtot, 0.0);
95 Array<OneD, NekDouble> out_dEta1 = out_dEta0 + Qtot;
96 Array<OneD, NekDouble> out_dEta2 = out_dEta1 + Qtot;
97
98 bool Do_2 = (out_dxi2.size() > 0) ? true : false;
99 bool Do_1 = (out_dxi1.size() > 0) ? true : false;
100
101 if (Do_2) // Need all local derivatives
102 {
103 PhysTensorDeriv(inarray, out_dEta0, out_dEta1, out_dEta2);
104 }
105 else if (Do_1) // Need 0 and 1 derivatives
106 {
107 PhysTensorDeriv(inarray, out_dEta0, out_dEta1, NullNekDouble1DArray);
108 }
109 else // Only need Eta0 derivaitve
110 {
111 PhysTensorDeriv(inarray, out_dEta0, NullNekDouble1DArray,
113 }
114
115 Array<OneD, const NekDouble> eta_0, eta_1, eta_2;
116 eta_0 = m_base[0]->GetZ();
117 eta_1 = m_base[1]->GetZ();
118 eta_2 = m_base[2]->GetZ();
119
120 // calculate 2.0/((1-eta_1)(1-eta_2)) Out_dEta0
121
122 NekDouble *dEta0 = &out_dEta0[0];
123 NekDouble fac;
124 for (int k = 0; k < Q2; ++k)
125 {
126 for (int j = 0; j < Q1; ++j, dEta0 += Q0)
127 {
128 Vmath::Smul(Q0, 2.0 / (1.0 - eta_1[j]), dEta0, 1, dEta0, 1);
129 }
130 fac = 1.0 / (1.0 - eta_2[k]);
131 Vmath::Smul(Q0 * Q1, fac, &out_dEta0[0] + k * Q0 * Q1, 1,
132 &out_dEta0[0] + k * Q0 * Q1, 1);
133 }
134
135 if (out_dxi0.size() > 0)
136 {
137 // out_dxi0 = 4.0/((1-eta_1)(1-eta_2)) Out_dEta0
138 Vmath::Smul(Qtot, 2.0, out_dEta0, 1, out_dxi0, 1);
139 }
140
141 if (Do_1 || Do_2)
142 {
143 Array<OneD, NekDouble> Fac0(Q0);
144 Vmath::Sadd(Q0, 1.0, eta_0, 1, Fac0, 1);
145
146 // calculate 2.0*(1+eta_0)/((1-eta_1)(1-eta_2)) Out_dEta0
147 for (int k = 0; k < Q1 * Q2; ++k)
148 {
149 Vmath::Vmul(Q0, &Fac0[0], 1, &out_dEta0[0] + k * Q0, 1,
150 &out_dEta0[0] + k * Q0, 1);
151 }
152 // calculate 2/(1.0-eta_2) out_dEta1
153 for (int k = 0; k < Q2; ++k)
154 {
155 Vmath::Smul(Q0 * Q1, 2.0 / (1.0 - eta_2[k]),
156 &out_dEta1[0] + k * Q0 * Q1, 1,
157 &out_dEta1[0] + k * Q0 * Q1, 1);
158 }
159
160 if (Do_1)
161 {
162 // calculate out_dxi1 = 2.0(1+eta_0)/((1-eta_1)(1-eta_2)) Out_dEta0
163 // + 2/(1.0-eta_2) out_dEta1
164 Vmath::Vadd(Qtot, out_dEta0, 1, out_dEta1, 1, out_dxi1, 1);
165 }
166
167 if (Do_2)
168 {
169 // calculate (1 + eta_1)/(1 -eta_2)*out_dEta1
170 NekDouble *dEta1 = &out_dEta1[0];
171 for (int k = 0; k < Q2; ++k)
172 {
173 for (int j = 0; j < Q1; ++j, dEta1 += Q0)
174 {
175 Vmath::Smul(Q0, (1.0 + eta_1[j]) / 2.0, dEta1, 1, dEta1, 1);
176 }
177 }
178
179 // calculate out_dxi2 =
180 // 2.0(1+eta_0)/((1-eta_1)(1-eta_2)) Out_dEta0 +
181 // (1 + eta_1)/(1 -eta_2)*out_dEta1 + out_dEta2
182 Vmath::Vadd(Qtot, out_dEta0, 1, out_dEta1, 1, out_dxi2, 1);
183 Vmath::Vadd(Qtot, out_dEta2, 1, out_dxi2, 1, out_dxi2, 1);
184 }
185 }
186}
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.
static Array< OneD, NekDouble > NullNekDouble1DArray
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition: Vmath.hpp:194

References Nektar::StdRegions::StdExpansion::m_base, Nektar::NullNekDouble1DArray, Nektar::StdRegions::StdExpansion3D::PhysTensorDeriv(), Vmath::Sadd(), Vmath::Smul(), Vmath::Vadd(), and Vmath::Vmul().

Referenced by v_PhysDeriv(), and v_StdPhysDeriv().

◆ v_PhysDeriv() [2/2]

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 194 of file StdTetExp.cpp.

197{
198 switch (dir)
199 {
200 case 0:
201 {
202 v_PhysDeriv(inarray, outarray, NullNekDouble1DArray,
204 break;
205 }
206 case 1:
207 {
208 v_PhysDeriv(inarray, NullNekDouble1DArray, outarray,
210 break;
211 }
212 case 2:
213 {
215 outarray);
216 break;
217 }
218 default:
219 {
220 ASSERTL1(false, "input dir is out of range");
221 }
222 break;
223 }
224}
void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dx, Array< OneD, NekDouble > &out_dy, Array< OneD, NekDouble > &out_dz) override
Calculate the derivative of the physical points.
Definition: StdTetExp.cpp:83

References ASSERTL1, Nektar::NullNekDouble1DArray, and v_PhysDeriv().

◆ v_PhysEvaluate()

NekDouble Nektar::StdRegions::StdTetExp::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::TetExp.

Definition at line 857 of file StdTetExp.cpp.

860{
861 // Collapse coordinates
862 Array<OneD, NekDouble> coll(3, 0.0);
863 LocCoordToLocCollapsed(coord, coll);
864
865 // If near singularity do the old interpolation matrix method
866 if ((1 - coll[1]) < 1e-5 || (1 - coll[2]) < 1e-5)
867 {
868 int totPoints = GetTotPoints();
869 Array<OneD, NekDouble> EphysDeriv0(totPoints), EphysDeriv1(totPoints),
870 EphysDeriv2(totPoints);
871 PhysDeriv(inarray, EphysDeriv0, EphysDeriv1, EphysDeriv2);
872
873 Array<OneD, DNekMatSharedPtr> I(3);
874 I[0] = GetBase()[0]->GetI(coll);
875 I[1] = GetBase()[1]->GetI(coll + 1);
876 I[2] = GetBase()[2]->GetI(coll + 2);
877
878 firstOrderDerivs[0] = PhysEvaluate(I, EphysDeriv0);
879 firstOrderDerivs[1] = PhysEvaluate(I, EphysDeriv1);
880 firstOrderDerivs[2] = PhysEvaluate(I, EphysDeriv2);
881 return PhysEvaluate(I, inarray);
882 }
883
884 std::array<NekDouble, 3> interDeriv;
885 NekDouble val = BaryTensorDeriv(coll, inarray, interDeriv);
886
887 // calculate 2.0/((1-eta_1)(1-eta_2)) * Out_dEta0
888 NekDouble temp = 2.0 / ((1 - coll[1]) * (1 - coll[2]));
889 interDeriv[0] *= temp;
890
891 // out_dxi0 = 4.0/((1-eta_1)(1-eta_2)) * Out_dEta0
892 firstOrderDerivs[0] = 2 * interDeriv[0];
893
894 // fac0 = 1 + eta_0
895 NekDouble fac0;
896 fac0 = 1 + coll[0];
897
898 // calculate 2.0*(1+eta_0)/((1-eta_1)(1-eta_2)) * Out_dEta0
899 interDeriv[0] *= fac0;
900
901 // calculate 2/(1.0-eta_2) * out_dEta1
902 fac0 = 2 / (1 - coll[2]);
903 interDeriv[1] *= fac0;
904
905 // calculate out_dxi1 = 2.0(1+eta_0)/((1-eta_1)(1-eta_2))
906 // * Out_dEta0 + 2/(1.0-eta_2) out_dEta1
907 firstOrderDerivs[1] = interDeriv[0] + interDeriv[1];
908
909 // calculate (1 + eta_1)/(1 -eta_2)*out_dEta1
910 fac0 = (1 + coll[1]) / 2;
911 interDeriv[1] *= fac0;
912
913 // calculate out_dxi2 =
914 // 2.0(1+eta_0)/((1-eta_1)(1-eta_2)) Out_dEta0 +
915 // (1 + eta_1)/(1 -eta_2)*out_dEta1 + out_dEta2
916 firstOrderDerivs[2] = interDeriv[0] + interDeriv[1] + interDeriv[2];
917
918 return val;
919}
NekDouble BaryTensorDeriv(const Array< OneD, NekDouble > &coord, const Array< OneD, const NekDouble > &inarray, std::array< NekDouble, 3 > &firstOrderDerivs)
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:134
const Array< OneD, const LibUtilities::BasisSharedPtr > & GetBase() const
This function gets the shared point to basis.
Definition: StdExpansion.h:100
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.
Definition: StdExpansion.h:919
void PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Definition: StdExpansion.h:849

References Nektar::StdRegions::StdExpansion3D::BaryTensorDeriv(), Nektar::StdRegions::StdExpansion::GetBase(), Nektar::StdRegions::StdExpansion::GetTotPoints(), Nektar::StdRegions::StdExpansion::LocCoordToLocCollapsed(), Nektar::StdRegions::StdExpansion::PhysDeriv(), and Nektar::StdRegions::StdExpansion::PhysEvaluate().

◆ v_PhysEvaluateBasis()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 815 of file StdTetExp.cpp.

817{
818 Array<OneD, NekDouble> coll(3);
819 LocCoordToLocCollapsed(coords, coll);
820
821 const int nm1 = m_base[1]->GetNumModes();
822 const int nm2 = m_base[2]->GetNumModes();
823
824 const int b = 2 * nm2 + 1;
825 const int mode0 = floor(0.5 * (b - sqrt(b * b - 8.0 * mode / nm1)));
826 const int tmp =
827 mode - nm1 * (mode0 * (nm2 - 1) + 1 - (mode0 - 2) * (mode0 - 1) / 2);
828 const int mode1 = tmp / (nm2 - mode0);
829 const int mode2 = tmp % (nm2 - mode0);
830
832 {
833 // Handle the collapsed vertices and edges in the modified
834 // basis.
835 if (mode == 1)
836 {
837 // Collapsed top vertex
838 return StdExpansion::BaryEvaluateBasis<2>(coll[2], 1);
839 }
840 else if (mode0 == 0 && mode2 == 1)
841 {
842 return StdExpansion::BaryEvaluateBasis<1>(coll[1], 0) *
843 StdExpansion::BaryEvaluateBasis<2>(coll[2], 1);
844 }
845 else if (mode0 == 1 && mode1 == 1 && mode2 == 0)
846 {
847 return StdExpansion::BaryEvaluateBasis<0>(coll[0], 0) *
848 StdExpansion::BaryEvaluateBasis<1>(coll[1], 1);
849 }
850 }
851
852 return StdExpansion::BaryEvaluateBasis<0>(coll[0], mode0) *
853 StdExpansion::BaryEvaluateBasis<1>(coll[1], mode1) *
854 StdExpansion::BaryEvaluateBasis<2>(coll[2], mode2);
855}
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::LocCoordToLocCollapsed(), Nektar::StdRegions::StdExpansion::m_base, and tinysimd::sqrt().

◆ v_ReduceOrderCoeffs()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2173 of file StdTetExp.cpp.

2176{
2177 int nquad0 = m_base[0]->GetNumPoints();
2178 int nquad1 = m_base[1]->GetNumPoints();
2179 int nquad2 = m_base[2]->GetNumPoints();
2180 int nqtot = nquad0 * nquad1 * nquad2;
2181 int nmodes0 = m_base[0]->GetNumModes();
2182 int nmodes1 = m_base[1]->GetNumModes();
2183 int nmodes2 = m_base[2]->GetNumModes();
2184 int numMax = nmodes0;
2185
2186 Array<OneD, NekDouble> coeff(m_ncoeffs);
2187 Array<OneD, NekDouble> coeff_tmp1(m_ncoeffs, 0.0);
2188 Array<OneD, NekDouble> coeff_tmp2(m_ncoeffs, 0.0);
2189 Array<OneD, NekDouble> phys_tmp(nqtot, 0.0);
2190 Array<OneD, NekDouble> tmp, tmp2, tmp3, tmp4;
2191
2192 Vmath::Vcopy(m_ncoeffs, inarray, 1, coeff_tmp2, 1);
2193
2194 const LibUtilities::PointsKey Pkey0 = m_base[0]->GetPointsKey();
2195 const LibUtilities::PointsKey Pkey1 = m_base[1]->GetPointsKey();
2196 const LibUtilities::PointsKey Pkey2 = m_base[2]->GetPointsKey();
2197
2198 LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A, nmodes0, Pkey0);
2199 LibUtilities::BasisKey bortho1(LibUtilities::eOrtho_B, nmodes1, Pkey1);
2200 LibUtilities::BasisKey bortho2(LibUtilities::eOrtho_C, nmodes2, Pkey2);
2201
2202 Vmath::Zero(m_ncoeffs, coeff_tmp2, 1);
2203
2206 bortho0, bortho1, bortho2);
2207
2208 BwdTrans(inarray, phys_tmp);
2209 OrthoTetExp->FwdTrans(phys_tmp, coeff);
2210
2211 Vmath::Zero(m_ncoeffs, outarray, 1);
2212
2213 // filtering
2214 int cnt = 0;
2215 for (int u = 0; u < numMin; ++u)
2216 {
2217 for (int i = 0; i < numMin - u; ++i)
2218 {
2219 Vmath::Vcopy(numMin - u - i, tmp = coeff + cnt, 1,
2220 tmp2 = coeff_tmp1 + cnt, 1);
2221 cnt += numMax - u - i;
2222 }
2223 for (int i = numMin; i < numMax - u; ++i)
2224 {
2225 cnt += numMax - u - i;
2226 }
2227 }
2228
2229 OrthoTetExp->BwdTrans(coeff_tmp1, phys_tmp);
2230 FwdTrans(phys_tmp, outarray);
2231}
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Backward transformation from coefficient space to physical space.
Definition: StdExpansion.h:424
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Forward transformation from physical space to coefficient space.
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:42
std::shared_ptr< StdTetExp > StdTetExpSharedPtr
Definition: StdTetExp.h:233
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::StdRegions::StdExpansion::BwdTrans(), Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::eOrtho_B, Nektar::LibUtilities::eOrtho_C, Nektar::StdRegions::StdExpansion::FwdTrans(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Vmath::Vcopy(), and Vmath::Zero().

◆ v_StdPhysDeriv() [1/2]

void Nektar::StdRegions::StdTetExp::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 226 of file StdTetExp.cpp.

230{
231 StdTetExp::v_PhysDeriv(inarray, out_d0, out_d1, out_d2);
232}

References v_PhysDeriv().

◆ v_StdPhysDeriv() [2/2]

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 234 of file StdTetExp.cpp.

237{
238 StdTetExp::v_PhysDeriv(dir, inarray, outarray);
239}

References v_PhysDeriv().

◆ v_SVVLaplacianFilter()

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TetExp.

Definition at line 2031 of file StdTetExp.cpp.

2033{
2034 // To do : 1) add a test to ensure 0 \leq SvvCutoff \leq 1.
2035 // 2) check if the transfer function needs an analytical
2036 // Fourier transform.
2037 // 3) if it doesn't : find a transfer function that renders
2038 // the if( cutoff_a ...) useless to reduce computational
2039 // cost.
2040 // 4) add SVVDiffCoef to both models!!
2041
2042 int qa = m_base[0]->GetNumPoints();
2043 int qb = m_base[1]->GetNumPoints();
2044 int qc = m_base[2]->GetNumPoints();
2045 int nmodes_a = m_base[0]->GetNumModes();
2046 int nmodes_b = m_base[1]->GetNumModes();
2047 int nmodes_c = m_base[2]->GetNumModes();
2048
2049 // Declare orthogonal basis.
2050 LibUtilities::PointsKey pa(qa, m_base[0]->GetPointsType());
2051 LibUtilities::PointsKey pb(qb, m_base[1]->GetPointsType());
2052 LibUtilities::PointsKey pc(qc, m_base[2]->GetPointsType());
2053
2054 LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A, nmodes_a, pa);
2055 LibUtilities::BasisKey Bb(LibUtilities::eOrtho_B, nmodes_b, pb);
2056 LibUtilities::BasisKey Bc(LibUtilities::eOrtho_C, nmodes_c, pc);
2057
2058 StdTetExp OrthoExp(Ba, Bb, Bc);
2059
2060 Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
2061 int i, j, k, cnt = 0;
2062
2063 // project onto physical space.
2064 OrthoExp.FwdTrans(array, orthocoeffs);
2065
2066 if (mkey.ConstFactorExists(eFactorSVVPowerKerDiffCoeff))
2067 {
2068 // Rodrigo's power kernel
2069 NekDouble cutoff = mkey.GetConstFactor(eFactorSVVCutoffRatio);
2070 NekDouble SvvDiffCoeff =
2071 mkey.GetConstFactor(eFactorSVVPowerKerDiffCoeff) *
2072 mkey.GetConstFactor(eFactorSVVDiffCoeff);
2073
2074 for (i = 0; i < nmodes_a; ++i)
2075 {
2076 for (j = 0; j < nmodes_b - j; ++j)
2077 {
2078 NekDouble fac1 = std::max(
2079 pow((1.0 * i) / (nmodes_a - 1), cutoff * nmodes_a),
2080 pow((1.0 * j) / (nmodes_b - 1), cutoff * nmodes_b));
2081
2082 for (k = 0; k < nmodes_c - i - j; ++k)
2083 {
2084 NekDouble fac =
2085 std::max(fac1, pow((1.0 * k) / (nmodes_c - 1),
2086 cutoff * nmodes_c));
2087
2088 orthocoeffs[cnt] *= SvvDiffCoeff * fac;
2089 cnt++;
2090 }
2091 }
2092 }
2093 }
2094 else if (mkey.ConstFactorExists(
2095 eFactorSVVDGKerDiffCoeff)) // Rodrigo/Mansoor's DG Kernel
2096 {
2097 NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDGKerDiffCoeff) *
2098 mkey.GetConstFactor(eFactorSVVDiffCoeff);
2099
2100 int max_abc = max(nmodes_a - kSVVDGFiltermodesmin,
2101 nmodes_b - kSVVDGFiltermodesmin);
2102 max_abc = max(max_abc, nmodes_c - kSVVDGFiltermodesmin);
2103 // clamp max_abc
2104 max_abc = max(max_abc, 0);
2105 max_abc = min(max_abc, kSVVDGFiltermodesmax - kSVVDGFiltermodesmin);
2106
2107 for (i = 0; i < nmodes_a; ++i)
2108 {
2109 for (j = 0; j < nmodes_b - j; ++j)
2110 {
2111 int maxij = max(i, j);
2112
2113 for (k = 0; k < nmodes_c - i - j; ++k)
2114 {
2115 int maxijk = max(maxij, k);
2116 maxijk = min(maxijk, kSVVDGFiltermodesmax - 1);
2117
2118 orthocoeffs[cnt] *=
2119 SvvDiffCoeff * kSVVDGFilter[max_abc][maxijk];
2120 cnt++;
2121 }
2122 }
2123 }
2124 }
2125 else
2126 {
2127
2128 // SVV filter paramaters (how much added diffusion
2129 // relative to physical one and fraction of modes from
2130 // which you start applying this added diffusion)
2131
2132 NekDouble SvvDiffCoeff =
2133 mkey.GetConstFactor(StdRegions::eFactorSVVDiffCoeff);
2134 NekDouble SVVCutOff =
2135 mkey.GetConstFactor(StdRegions::eFactorSVVCutoffRatio);
2136
2137 // Defining the cut of mode
2138 int cutoff_a = (int)(SVVCutOff * nmodes_a);
2139 int cutoff_b = (int)(SVVCutOff * nmodes_b);
2140 int cutoff_c = (int)(SVVCutOff * nmodes_c);
2141 int nmodes = min(min(nmodes_a, nmodes_b), nmodes_c);
2142 NekDouble cutoff = min(min(cutoff_a, cutoff_b), cutoff_c);
2143 NekDouble epsilon = 1;
2144
2145 //------"New" Version August 22nd '13--------------------
2146 for (i = 0; i < nmodes_a; ++i)
2147 {
2148 for (j = 0; j < nmodes_b - i; ++j)
2149 {
2150 for (k = 0; k < nmodes_c - i - j; ++k)
2151 {
2152 if (i + j + k >= cutoff)
2153 {
2154 orthocoeffs[cnt] *= ((SvvDiffCoeff)*exp(
2155 -(i + j + k - nmodes) * (i + j + k - nmodes) /
2156 ((NekDouble)((i + j + k - cutoff + epsilon) *
2157 (i + j + k - cutoff + epsilon)))));
2158 }
2159 else
2160 {
2161 orthocoeffs[cnt] *= 0.0;
2162 }
2163 cnt++;
2164 }
2165 }
2166 }
2167 }
2168
2169 // backward transform to physical space
2170 OrthoExp.BwdTrans(orthocoeffs, array);
2171}
StdTetExp(const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
Definition: StdTetExp.cpp:42
const int kSVVDGFiltermodesmin
Definition: StdRegions.hpp:500
const int kSVVDGFiltermodesmax
Definition: StdRegions.hpp:501
const NekDouble kSVVDGFilter[9][11]
Definition: StdRegions.hpp:503

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::LibUtilities::eOrtho_B, Nektar::LibUtilities::eOrtho_C, 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.