Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Protected Member Functions | List of all members
Nektar::StdRegions::StdHexExp Class Reference

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

#include <StdHexExp.h>

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

Public Member Functions

 StdHexExp ()
 
 StdHexExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
 
 StdHexExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc, NekDouble *coeffs, NekDouble *phys)
 
 StdHexExp (const StdHexExp &T)
 
 ~StdHexExp ()
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion3D
 StdExpansion3D ()
 
 StdExpansion3D (int numcoeffs, const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
 
 StdExpansion3D (const StdExpansion3D &T)
 
virtual ~StdExpansion3D ()
 
void PhysTensorDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d1, Array< OneD, NekDouble > &outarray_d2, Array< OneD, NekDouble > &outarray_d3)
 Calculate the 3D derivative in the local tensor/collapsed coordinate at the physical points. More...
 
void BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
 
void IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
 
- 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::BasisSharedPtr
GetBasis (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 GetNedges () const
 This function returns the number of edges of the expansion domain. More...
 
int GetEdgeNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th edge. More...
 
int GetTotalEdgeIntNcoeffs () const
 
int GetEdgeNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th edge. More...
 
int DetCartesianDirOfEdge (const int edge)
 
const LibUtilities::BasisKey DetEdgeBasisKey (const int i) const
 
const LibUtilities::BasisKey DetFaceBasisKey (const int i, const int k) const
 
int GetFaceNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th face. More...
 
int GetFaceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th face. More...
 
int GetFaceIntNcoeffs (const int i) const
 
int GetTotalFaceIntNcoeffs () const
 
int GetTraceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th edge/face. More...
 
LibUtilities::PointsKey GetFacePointsKey (const int i, const int j) const
 
int NumBndryCoeffs (void) const
 
int NumDGBndryCoeffs (void) const
 
LibUtilities::BasisType GetEdgeBasisType (const int i) const
 This function returns the type of expansion basis on the i-th edge. More...
 
const LibUtilities::PointsKey GetNodalPointsKey () const
 This function returns the type of expansion Nodal point type if defined. More...
 
int GetNfaces () const
 This function returns the number of faces of the expansion domain. More...
 
int GetNtrace () 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...
 
boost::shared_ptr< StdExpansionGetStdExp (void) const
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion ()
 
bool IsNodalNonTensorialExp ()
 
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space. More...
 
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space. More...
 
void FwdTrans_BndConstrained (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)
 
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)
 
IndexMapValuesSharedPtr GetIndexMap (const IndexMapKey &ikey)
 
const Array< OneD, const
NekDouble > & 
GetPhysNormals (void)
 
void SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
virtual void SetUpPhysNormals (const int edge)
 
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)
 
StdRegions::Orientation GetForient (int face)
 
StdRegions::Orientation GetEorient (int edge)
 
StdRegions::Orientation GetPorient (int point)
 
StdRegions::Orientation GetCartesianEorient (int edge)
 
void SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
 
void SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
int CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
void ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs)
 
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 GetEdgeInteriorMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
 
void GetFaceInteriorMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
 
void GetEdgeToElementMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int P=-1)
 
void GetFaceToElementMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int nummodesA=-1, int nummodesB=-1)
 
void GetEdgePhysVals (const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Extract the physical values along edge edge from inarray into outarray following the local edge orientation and point distribution defined by defined in EdgeExp. More...
 
void GetEdgePhysVals (const int edge, const boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetTracePhysVals (const int edge, const boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetVertexPhysVals (const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray)
 
void GetEdgeInterpVals (const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetEdgeQFactors (const int edge, Array< OneD, NekDouble > &outarray)
 Extract the metric factors to compute the contravariant fluxes along edge edge and stores them into outarray following the local edge orientation (i.e. anticlockwise convention). More...
 
void GetFacePhysVals (const int face, const boost::shared_ptr< StdExpansion > &FaceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=eNoOrientation)
 
void GetEdgePhysMap (const int edge, Array< OneD, int > &outarray)
 
void GetFacePhysMap (const int face, Array< OneD, int > &outarray)
 
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 LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
DNekMatSharedPtr GenMatrix (const StdMatrixKey &mkey)
 
void PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void PhysDeriv_s (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds)
 
void PhysDeriv_n (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn)
 
void PhysDirectionalDeriv (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &outarray)
 
void StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void AddRobinMassMatrix (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
void AddRobinEdgeContribution (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, Array< OneD, NekDouble > &coeffs)
 
NekDouble PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta. More...
 
const boost::shared_ptr
< SpatialDomains::GeomFactors > & 
GetMetricInfo (void) const
 
virtual int v_GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
virtual const Array< OneD,
const NekDouble > & 
v_GetPhysNormals (void)
 
virtual void v_SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
virtual void v_SetUpPhysNormals (const int edge)
 
virtual void v_ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmode_offset, NekDouble *coeffs)
 Unpack data from input file assuming it comes from the same expansion type. More...
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
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)
 
virtual StdRegions::Orientation v_GetForient (int face)
 
virtual StdRegions::Orientation v_GetEorient (int edge)
 
virtual StdRegions::Orientation v_GetCartesianEorient (int edge)
 
virtual StdRegions::Orientation v_GetPorient (int point)
 
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 NormalVectorGetEdgeNormal (const int edge) const
 
void ComputeEdgeNormal (const int edge)
 
void NegateEdgeNormal (const int edge)
 
bool EdgeNormalNegated (const int edge)
 
void ComputeFaceNormal (const int face)
 
void NegateFaceNormal (const int face)
 
bool FaceNormalNegated (const int face)
 
void ComputeVertexNormal (const int vertex)
 
const NormalVectorGetFaceNormal (const int face) const
 
const NormalVectorGetVertexNormal (const int vertex) const
 
const NormalVectorGetSurfaceNormal (const int id) const
 
const LibUtilities::PointsKeyVector GetPointsKeys () const
 
Array< OneD, unsigned int > GetEdgeInverseBoundaryMap (int eid)
 
Array< OneD, unsigned int > GetFaceInverseBoundaryMap (int fid, StdRegions::Orientation faceOrient=eNoOrientation)
 
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 >
boost::shared_ptr< T > as ()
 
void IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 

Protected Member Functions

virtual void v_PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
 Differentiation Methods. More...
 
virtual void v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
 
virtual void v_StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
 
virtual void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTBase_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multbyweights=true)
 
virtual void v_IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
 
virtual void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTDerivBase_MatOp (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 
virtual void v_FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 
virtual int v_GetNverts () const
 
virtual int v_GetNedges () const
 
virtual int v_GetNfaces () const
 
virtual LibUtilities::ShapeType v_DetShapeType () const
 
virtual int v_NumBndryCoeffs () const
 
virtual int v_NumDGBndryCoeffs () const
 
virtual int v_GetEdgeNcoeffs (const int i) const
 
virtual int v_GetTotalEdgeIntNcoeffs () const
 
virtual int v_GetFaceNcoeffs (const int i) const
 
virtual int v_GetFaceIntNcoeffs (const int i) const
 
virtual int v_GetTotalFaceIntNcoeffs () const
 
virtual int v_GetFaceNumPoints (const int i) const
 
virtual LibUtilities::PointsKey v_GetFacePointsKey (const int i, const int j) const
 
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
virtual const
LibUtilities::BasisKey 
v_DetFaceBasisKey (const int i, const int k) const
 
virtual LibUtilities::BasisType v_GetEdgeBasisType (const int i) const
 
virtual bool v_IsBoundaryInteriorExpansion ()
 
virtual void v_GetCoords (Array< OneD, NekDouble > &coords_x, Array< OneD, NekDouble > &coords_y, Array< OneD, NekDouble > &coords_z)
 
virtual void v_GetFaceToElementMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int nummodesA=-1, int nummodesB=-1)
 
virtual int v_GetVertexMap (int localVertexId, bool useCoeffPacking=false)
 
virtual void v_GetEdgeInteriorMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
 
virtual void v_GetFaceInteriorMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
 
virtual void v_GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
virtual void v_GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
virtual DNekMatSharedPtr v_GenMatrix (const StdMatrixKey &mkey)
 
virtual DNekMatSharedPtr v_CreateStdMatrix (const StdMatrixKey &mkey)
 
virtual void v_MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_GeneralMatrixOp_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion3D
virtual NekDouble v_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...
 
virtual NekDouble v_PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 
virtual void v_LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
 
virtual void v_HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
 
virtual NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray)
 Integrates the specified function over the domain. More...
 
virtual void v_NegateFaceNormal (const int face)
 
virtual bool v_FaceNormalNegated (const int face)
 
virtual int v_GetTraceNcoeffs (const int i) const
 
- 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...
 
IndexMapValuesSharedPtr CreateIndexMap (const IndexMapKey &ikey)
 Create an IndexMap which contains mapping information linking any specific element shape with either its boundaries, edges, faces, verteces, etc. More...
 
void BwdTrans_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GeneralMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
 
void LaplacianMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp_MatFree (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void HelmholtzMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_SetCoeffsToOrientation (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)
 

Additional Inherited Members

- Protected Attributes inherited from Nektar::StdRegions::StdExpansion3D
std::map< int, NormalVectorm_faceNormals
 
std::map< int, bool > m_negatedNormals
 
- Protected Attributes inherited from Nektar::StdRegions::StdExpansion
Array< OneD,
LibUtilities::BasisSharedPtr
m_base
 
int m_elmt_id
 
int m_ncoeffs
 
LibUtilities::NekManager
< StdMatrixKey, DNekMat,
StdMatrixKey::opLess
m_stdMatrixManager
 
LibUtilities::NekManager
< StdMatrixKey, DNekBlkMat,
StdMatrixKey::opLess
m_stdStaticCondMatrixManager
 
LibUtilities::NekManager
< IndexMapKey, IndexMapValues,
IndexMapKey::opLess
m_IndexMapManager
 

Detailed Description

Class representing a hexehedral element in reference space.

Definition at line 48 of file StdHexExp.h.

Constructor & Destructor Documentation

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

Definition at line 49 of file StdHexExp.cpp.

50  {
51  }
Nektar::StdRegions::StdHexExp::StdHexExp ( const LibUtilities::BasisKey Ba,
const LibUtilities::BasisKey Bb,
const LibUtilities::BasisKey Bc 
)

Definition at line 54 of file StdHexExp.cpp.

56  :
57  StdExpansion(Ba.GetNumModes()*Bb.GetNumModes()*Bc.GetNumModes(), 3,
58  Ba, Bb, Bc),
59  StdExpansion3D(Ba.GetNumModes()*Bb.GetNumModes()*Bc.GetNumModes(),
60  Ba, Bb, Bc)
61  {
62  }
StdExpansion()
Default Constructor.
Nektar::StdRegions::StdHexExp::StdHexExp ( const LibUtilities::BasisKey Ba,
const LibUtilities::BasisKey Bb,
const LibUtilities::BasisKey Bc,
NekDouble coeffs,
NekDouble phys 
)

Definition at line 65 of file StdHexExp.cpp.

70  {
71  }
Nektar::StdRegions::StdHexExp::StdHexExp ( const StdHexExp T)

Definition at line 74 of file StdHexExp.cpp.

74  :
75  StdExpansion(T),
77  {
78  }
StdExpansion()
Default Constructor.
Nektar::StdRegions::StdHexExp::~StdHexExp ( )

Definition at line 81 of file StdHexExp.cpp.

82  {
83  }

Member Function Documentation

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

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

Parameters
inarray?
outarray?

Implements Nektar::StdRegions::StdExpansion.

Definition at line 184 of file StdHexExp.cpp.

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

187  {
190  "Basis[1] is not a general tensor type");
191 
194  "Basis[2] is not a general tensor type");
195 
196  if(m_base[0]->Collocation() && m_base[1]->Collocation()
197  && m_base[2]->Collocation())
198  {
200  * m_base[1]->GetNumPoints()
201  * m_base[2]->GetNumPoints(),
202  inarray, 1, outarray, 1);
203  }
204  else
205  {
206  StdHexExp::BwdTrans_SumFac(inarray,outarray);
207  }
208  }
Principle Modified Functions .
Definition: BasisType.h:51
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:229
Principle Orthogonal Functions .
Definition: BasisType.h:47
Principle Modified Functions .
Definition: BasisType.h:50
Principle Orthogonal Functions .
Definition: BasisType.h:48
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
void BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void Nektar::StdRegions::StdHexExp::v_BwdTrans_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 214 of file StdHexExp.cpp.

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

216  {
217  Array<OneD, NekDouble> wsp(m_base[0]->GetNumPoints()*
218  m_base[2]->GetNumModes()*
219  (m_base[1]->GetNumModes() + m_base[1]->GetNumPoints())); // FIX THIS
220 
221  BwdTrans_SumFacKernel(m_base[0]->GetBdata(),
222  m_base[1]->GetBdata(),
223  m_base[2]->GetBdata(),
224  inarray,outarray,wsp,true,true,true);
225  }
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)
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:229
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdHexExp::v_BwdTrans_SumFacKernel ( const Array< OneD, const NekDouble > &  base0,
const Array< OneD, const NekDouble > &  base1,
const Array< OneD, const NekDouble > &  base2,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
Array< OneD, NekDouble > &  wsp,
bool  doCheckCollDir0,
bool  doCheckCollDir1,
bool  doCheckCollDir2 
)
protectedvirtual
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 241 of file StdHexExp.cpp.

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

251  {
252  int nquad0 = m_base[0]->GetNumPoints();
253  int nquad1 = m_base[1]->GetNumPoints();
254  int nquad2 = m_base[2]->GetNumPoints();
255  int nmodes0 = m_base[0]->GetNumModes();
256  int nmodes1 = m_base[1]->GetNumModes();
257  int nmodes2 = m_base[2]->GetNumModes();
258 
259  // Check if using collocation, if requested.
260  bool colldir0 = doCheckCollDir0?(m_base[0]->Collocation()):false;
261  bool colldir1 = doCheckCollDir1?(m_base[1]->Collocation()):false;
262  bool colldir2 = doCheckCollDir2?(m_base[2]->Collocation()):false;
263 
264  // If collocation in all directions, Physical values at quadrature
265  // points is just a copy of the modes.
266  if(colldir0 && colldir1 && colldir2)
267  {
268  Vmath::Vcopy(m_ncoeffs,inarray.get(),1,outarray.get(),1);
269  }
270  else
271  {
272  // Check sufficiently large workspace.
273  ASSERTL1(wsp.num_elements()>=nquad0*nmodes2*(nmodes1+nquad1),
274  "Workspace size is not sufficient");
275 
276  // Assign second half of workspace for 2nd DGEMM operation.
277  Array<OneD, NekDouble> wsp2 = wsp + nquad0*nmodes1*nmodes2;
278 
279  // BwdTrans in each direction using DGEMM
280  Blas::Dgemm('T','T', nmodes1*nmodes2, nquad0, nmodes0,
281  1.0, &inarray[0], nmodes0,
282  base0.get(), nquad0,
283  0.0, &wsp[0], nmodes1*nmodes2);
284  Blas::Dgemm('T','T', nquad0*nmodes2, nquad1, nmodes1,
285  1.0, &wsp[0], nmodes1,
286  base1.get(), nquad1,
287  0.0, &wsp2[0], nquad0*nmodes2);
288  Blas::Dgemm('T','T', nquad0*nquad1, nquad2, nmodes2,
289  1.0, &wsp2[0], nmodes2,
290  base2.get(), nquad2,
291  0.0, &outarray[0], nquad0*nquad1);
292  }
293  }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
int Nektar::StdRegions::StdHexExp::v_CalcNumberOfCoefficients ( const std::vector< unsigned int > &  nummodes,
int &  modes_offset 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 848 of file StdHexExp.cpp.

849  {
850  int nmodes = nummodes[modes_offset]*nummodes[modes_offset+1]*nummodes[modes_offset+2];
851  modes_offset += 3;
852 
853  return nmodes;
854  }
DNekMatSharedPtr Nektar::StdRegions::StdHexExp::v_CreateStdMatrix ( const StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 2219 of file StdHexExp.cpp.

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

2220  {
2221  return StdExpansion::CreateGeneralMatrix(mkey);
2222  }
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
const LibUtilities::BasisKey Nektar::StdRegions::StdHexExp::v_DetFaceBasisKey ( const int  i,
const int  k 
) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 857 of file StdHexExp.cpp.

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

859  {
860  ASSERTL2(i >= 0 && i <= 5, "face id is out of range");
861  ASSERTL2(k >= 0 && k <= 1, "basis key id is out of range");
862 
863  int dir = k;
864  switch(i)
865  {
866  case 0:
867  case 5:
868  dir = k;
869  break;
870  case 1:
871  case 3:
872  dir = 2*k;
873  break;
874  case 2:
875  case 4:
876  dir = k+1;
877  break;
878  }
879 
880  return EvaluateQuadFaceBasisKey(k,
881  m_base[dir]->GetBasisType(),
882  m_base[dir]->GetNumPoints(),
883  m_base[dir]->GetNumModes());
884  }
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:229
LibUtilities::BasisKey EvaluateQuadFaceBasisKey(const int facedir, const LibUtilities::BasisType faceDirBasisType, const int numpoints, const int nummodes)
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:240
Array< OneD, LibUtilities::BasisSharedPtr > m_base
LibUtilities::ShapeType Nektar::StdRegions::StdHexExp::v_DetShapeType ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 693 of file StdHexExp.cpp.

References Nektar::LibUtilities::eHexahedron.

void Nektar::StdRegions::StdHexExp::v_FillMode ( const int  mode,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual
Note
for hexahedral expansions _base0 modes run fastest.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 624 of file StdHexExp.cpp.

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

626  {
627  int i,j;
628  int nquad0 = m_base[0]->GetNumPoints();
629  int nquad1 = m_base[1]->GetNumPoints();
630  int nquad2 = m_base[2]->GetNumPoints();
631 
632  Array<OneD, const NekDouble> base0 = m_base[0]->GetBdata();
633  Array<OneD, const NekDouble> base1 = m_base[1]->GetBdata();
634  Array<OneD, const NekDouble> base2 = m_base[2]->GetBdata();
635 
636  int btmp0 = m_base[0]->GetNumModes();
637  int btmp1 = m_base[1]->GetNumModes();
638  int mode2 = mode/(btmp0*btmp1);
639  int mode1 = (mode-mode2*btmp0*btmp1)/btmp0;
640  int mode0 = (mode-mode2*btmp0*btmp1)%btmp0;
641 
642  ASSERTL2(mode2 == (int)floor((1.0*mode)/(btmp0*btmp1)),
643  "Integer Truncation not Equiv to Floor");
644  ASSERTL2(mode1 == (int)floor((1.0*mode-mode2*btmp0*btmp1)
645  /(btmp0*btmp1)),
646  "Integer Truncation not Equiv to Floor");
647  ASSERTL2(m_ncoeffs <= mode,
648  "calling argument mode is larger than total expansion "
649  "order");
650 
651  for(i = 0; i < nquad1*nquad2; ++i)
652  {
653  Vmath::Vcopy(nquad0,(NekDouble *)(base0.get() + mode0*nquad0),1,
654  &outarray[0]+i*nquad0, 1);
655  }
656 
657  for(j = 0; j < nquad2; ++j)
658  {
659  for(i = 0; i < nquad0; ++i)
660  {
661  Vmath::Vmul(nquad1,(NekDouble *)(base1.get() + mode1*nquad1),1,
662  &outarray[0]+i+j*nquad0*nquad1, nquad0,
663  &outarray[0]+i+j*nquad0*nquad1, nquad0);
664  }
665  }
666 
667  for(i = 0; i < nquad2; i++)
668  {
669  Blas::Dscal(nquad0*nquad1,base2[mode2*nquad2+i],
670  &outarray[0]+i*nquad0*nquad1,1);
671  }
672  }
double NekDouble
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:240
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:169
void Nektar::StdRegions::StdHexExp::v_FwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

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

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

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 305 of file StdHexExp.cpp.

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

308  {
309  // If using collocation expansion, coefficients match physical
310  // data points so just do a direct copy.
311  if( (m_base[0]->Collocation())
312  &&(m_base[1]->Collocation())
313  &&(m_base[2]->Collocation()) )
314  {
315  Vmath::Vcopy(GetNcoeffs(), &inarray[0], 1, &outarray[0], 1);
316  }
317  else
318  {
319  // Compute B^TWu
320  IProductWRTBase(inarray,outarray);
321 
322  // get Mass matrix inverse
323  StdMatrixKey masskey(eInvMass,DetShapeType(),*this);
324  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
325 
326  // copy inarray in case inarray == outarray
327  DNekVec in (m_ncoeffs,outarray);
328  DNekVec out(m_ncoeffs,outarray,eWrapper);
329 
330  // Solve for coefficients.
331  out = (*matsys)*in;
332 
333  }
334  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:470
void IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
this function calculates the inner product of a given function f with the different modes of the expa...
Definition: StdExpansion.h:629
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:700
NekVector< NekDouble > DNekVec
Definition: NekTypeDefs.hpp:49
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:131
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
void Nektar::StdRegions::StdHexExp::v_GeneralMatrixOp_MatOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdMatrixKey mkey 
)
protectedvirtual

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 2271 of file StdHexExp.cpp.

References Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::StdRegions::StdExpansion::m_stdMatrixManager, and Vmath::Vcopy().

2275  {
2277 
2278  if(inarray.get() == outarray.get())
2279  {
2280  Array<OneD,NekDouble> tmp(m_ncoeffs);
2281  Vmath::Vcopy(m_ncoeffs,inarray.get(),1,tmp.get(),1);
2282 
2283  Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, 1.0, mat->GetPtr().get(),
2284  m_ncoeffs, tmp.get(), 1, 0.0, outarray.get(), 1);
2285  }
2286  else
2287  {
2288  Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, 1.0, mat->GetPtr().get(),
2289  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
2290  }
2291  }
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
LibUtilities::NekManager< StdMatrixKey, DNekMat, StdMatrixKey::opLess > m_stdMatrixManager
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
DNekMatSharedPtr Nektar::StdRegions::StdHexExp::v_GenMatrix ( const StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 2213 of file StdHexExp.cpp.

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

2214  {
2215  return StdExpansion::CreateGeneralMatrix(mkey);
2216  }
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
void Nektar::StdRegions::StdHexExp::v_GetBoundaryMap ( Array< OneD, unsigned int > &  outarray)
protectedvirtual
Parameters
outarrayStorage for computed map.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2122 of file StdHexExp.cpp.

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

2123  {
2126  "BasisType is not a boundary interior form");
2129  "BasisType is not a boundary interior form");
2132  "BasisType is not a boundary interior form");
2133 
2134  int i;
2135  int nummodes [3] = {m_base[0]->GetNumModes(),
2136  m_base[1]->GetNumModes(),
2137  m_base[2]->GetNumModes()};
2138 
2139  int nBndCoeffs = NumBndryCoeffs();
2140 
2141  if(outarray.num_elements()!=nBndCoeffs)
2142  {
2143  outarray = Array<OneD, unsigned int>(nBndCoeffs);
2144  }
2145 
2146  const LibUtilities::BasisType Btype [3] = {GetBasisType(0),
2147  GetBasisType(1),
2148  GetBasisType(2)};
2149 
2150  int p,q,r;
2151  int cnt = 0;
2152 
2153  int BndIdx [3][2];
2154  int IntIdx [3][2];
2155 
2156  for(i = 0; i < 3; i++)
2157  {
2158  BndIdx[i][0] = 0;
2159 
2160  if( Btype[i] == LibUtilities::eModified_A)
2161  {
2162  BndIdx[i][1] = 1;
2163  IntIdx[i][0] = 2;
2164  IntIdx[i][1] = nummodes[i];
2165  }
2166  else
2167  {
2168  BndIdx[i][1] = nummodes[i]-1;
2169  IntIdx[i][0] = 1;
2170  IntIdx[i][1] = nummodes[i]-1;
2171  }
2172  }
2173 
2174 
2175  for(i = 0; i < 2; i++)
2176  {
2177  r = BndIdx[2][i];
2178  for( q = 0; q < nummodes[1]; q++)
2179  {
2180  for( p = 0; p < nummodes[0]; p++)
2181  {
2182  outarray[cnt++] = r*nummodes[0]*nummodes[1]+q*nummodes[0] + p;
2183  }
2184  }
2185  }
2186 
2187  for(r = IntIdx[2][0]; r < IntIdx[2][1]; r++)
2188  {
2189  for( i = 0; i < 2; i++)
2190  {
2191  q = BndIdx[1][i];
2192  for( p = 0; p < nummodes[0]; p++)
2193  {
2194  outarray[cnt++] = r*nummodes[0]*nummodes[1] +
2195  q*nummodes[0] + p;
2196  }
2197  }
2198 
2199  for( q = IntIdx[1][0]; q < IntIdx[1][1]; q++)
2200  {
2201  for( i = 0; i < 2; i++)
2202  {
2203  p = BndIdx[0][i];
2204  outarray[cnt++] = r*nummodes[0]*nummodes[1] +
2205  q*nummodes[0] + p;
2206  }
2207  }
2208  }
2209 
2210  sort(outarray.get(), outarray.get() + nBndCoeffs);
2211  }
Principle Modified Functions .
Definition: BasisType.h:49
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
Lagrange for SEM basis .
Definition: BasisType.h:53
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdHexExp::v_GetCoords ( Array< OneD, NekDouble > &  coords_x,
Array< OneD, NekDouble > &  coords_y,
Array< OneD, NekDouble > &  coords_z 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 904 of file StdHexExp.cpp.

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

907  {
908  Array<OneD, const NekDouble> eta_x = m_base[0]->GetZ();
909  Array<OneD, const NekDouble> eta_y = m_base[1]->GetZ();
910  Array<OneD, const NekDouble> eta_z = m_base[2]->GetZ();
911  int Qx = GetNumPoints(0);
912  int Qy = GetNumPoints(1);
913  int Qz = GetNumPoints(2);
914 
915  // Convert collapsed coordinates into cartesian coordinates:
916  // eta --> xi
917  for( int k = 0; k < Qz; ++k ) {
918  for( int j = 0; j < Qy; ++j ) {
919  for( int i = 0; i < Qx; ++i ) {
920  int s = i + Qx*(j + Qy*k);
921  xi_x[s] = eta_x[i];
922  xi_y[s] = eta_y[j];
923  xi_z[s] = eta_z[k];
924 
925  }
926  }
927  }
928  }
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:229
Array< OneD, LibUtilities::BasisSharedPtr > m_base
LibUtilities::BasisType Nektar::StdRegions::StdHexExp::v_GetEdgeBasisType ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 886 of file StdHexExp.cpp.

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

887  {
888  ASSERTL2((i >= 0)&&(i <= 11),"edge id is out of range");
889 
890  if((i == 0)||(i == 2)||(i==8)||(i==10))
891  {
892  return GetBasisType(0);
893  }
894  else if((i == 1)||(i == 3)||(i == 9)||(i == 11))
895  {
896  return GetBasisType(1);
897  }
898  else
899  {
900  return GetBasisType(2);
901  }
902  }
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:240
void Nektar::StdRegions::StdHexExp::v_GetEdgeInteriorMap ( const int  eid,
const Orientation  edgeOrient,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray 
)
protectedvirtual
Parameters
eidThe edge to compute the numbering for.
edgeOrientOrientation of the edge.
maparrayStorage for computed mapping array.
signarray?

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1420 of file StdHexExp.cpp.

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

1424  {
1427  "BasisType is not a boundary interior form");
1430  "BasisType is not a boundary interior form");
1433  "BasisType is not a boundary interior form");
1434 
1435  ASSERTL1((eid>=0)&&(eid<12),
1436  "local edge id must be between 0 and 11");
1437 
1438  int nEdgeIntCoeffs = GetEdgeNcoeffs(eid)-2;
1439 
1440  if(maparray.num_elements()!=nEdgeIntCoeffs)
1441  {
1442  maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
1443  }
1444 
1445  if(signarray.num_elements() != nEdgeIntCoeffs)
1446  {
1447  signarray = Array<OneD, int>(nEdgeIntCoeffs,1);
1448  }
1449  else
1450  {
1451  fill( signarray.get() , signarray.get()+nEdgeIntCoeffs, 1 );
1452  }
1453 
1454  int nummodes [3] = {m_base[0]->GetNumModes(),
1455  m_base[1]->GetNumModes(),
1456  m_base[2]->GetNumModes()};
1457 
1458  const LibUtilities::BasisType bType [3] = {GetBasisType(0),
1459  GetBasisType(1),
1460  GetBasisType(2)};
1461 
1462  bool reverseOrdering = false;
1463  bool signChange = false;
1464 
1465  int IdxRange [3][2] = {{0,0},{0,0},{0,0}};
1466 
1467  switch(eid)
1468  {
1469  case 0:
1470  case 1:
1471  case 2:
1472  case 3:
1473  {
1474  IdxRange[2][0] = 0;
1475  IdxRange[2][1] = 1;
1476  }
1477  break;
1478  case 8:
1479  case 9:
1480  case 10:
1481  case 11:
1482  {
1483  if( bType[2] == LibUtilities::eGLL_Lagrange)
1484  {
1485  IdxRange[2][0] = nummodes[2] - 1;
1486  IdxRange[2][1] = nummodes[2];
1487  }
1488  else
1489  {
1490  IdxRange[2][0] = 1;
1491  IdxRange[2][1] = 2;
1492  }
1493  }
1494  break;
1495  case 4:
1496  case 5:
1497  case 6:
1498  case 7:
1499  {
1500  if( bType[2] == LibUtilities::eGLL_Lagrange)
1501  {
1502  IdxRange[2][0] = 1;
1503  IdxRange[2][1] = nummodes[2] - 1;
1504 
1505  if(edgeOrient==eBackwards)
1506  {
1507  reverseOrdering = true;
1508  }
1509  }
1510  else
1511  {
1512  IdxRange[2][0] = 2;
1513  IdxRange[2][1] = nummodes[2];
1514 
1515  if(edgeOrient==eBackwards)
1516  {
1517  signChange = true;
1518  }
1519  }
1520  }
1521  break;
1522  }
1523 
1524  switch(eid)
1525  {
1526  case 0:
1527  case 4:
1528  case 5:
1529  case 8:
1530  {
1531  IdxRange[1][0] = 0;
1532  IdxRange[1][1] = 1;
1533  }
1534  break;
1535  case 2:
1536  case 6:
1537  case 7:
1538  case 10:
1539  {
1540  if( bType[1] == LibUtilities::eGLL_Lagrange)
1541  {
1542  IdxRange[1][0] = nummodes[1] - 1;
1543  IdxRange[1][1] = nummodes[1];
1544  }
1545  else
1546  {
1547  IdxRange[1][0] = 1;
1548  IdxRange[1][1] = 2;
1549  }
1550  }
1551  break;
1552  case 1:
1553  case 9:
1554  {
1555  if( bType[1] == LibUtilities::eGLL_Lagrange)
1556  {
1557  IdxRange[1][0] = 1;
1558  IdxRange[1][1] = nummodes[1] - 1;
1559 
1560  if(edgeOrient==eBackwards)
1561  {
1562  reverseOrdering = true;
1563  }
1564  }
1565  else
1566  {
1567  IdxRange[1][0] = 2;
1568  IdxRange[1][1] = nummodes[1];
1569 
1570  if(edgeOrient==eBackwards)
1571  {
1572  signChange = true;
1573  }
1574  }
1575  }
1576  break;
1577  case 3:
1578  case 11:
1579  {
1580  if( bType[1] == LibUtilities::eGLL_Lagrange)
1581  {
1582  IdxRange[1][0] = 1;
1583  IdxRange[1][1] = nummodes[1] - 1;
1584 
1585  if(edgeOrient==eForwards)
1586  {
1587  reverseOrdering = true;
1588  }
1589  }
1590  else
1591  {
1592  IdxRange[1][0] = 2;
1593  IdxRange[1][1] = nummodes[1];
1594 
1595  if(edgeOrient==eForwards)
1596  {
1597  signChange = true;
1598  }
1599  }
1600  }
1601  break;
1602  }
1603 
1604  switch(eid)
1605  {
1606  case 3:
1607  case 4:
1608  case 7:
1609  case 11:
1610  {
1611  IdxRange[0][0] = 0;
1612  IdxRange[0][1] = 1;
1613  }
1614  break;
1615  case 1:
1616  case 5:
1617  case 6:
1618  case 9:
1619  {
1620  if( bType[0] == LibUtilities::eGLL_Lagrange)
1621  {
1622  IdxRange[0][0] = nummodes[0] - 1;
1623  IdxRange[0][1] = nummodes[0];
1624  }
1625  else
1626  {
1627  IdxRange[0][0] = 1;
1628  IdxRange[0][1] = 2;
1629  }
1630  }
1631  break;
1632  case 0:
1633  case 8:
1634  {
1635  if( bType[0] == LibUtilities::eGLL_Lagrange)
1636  {
1637  IdxRange[0][0] = 1;
1638  IdxRange[0][1] = nummodes[0] - 1;
1639 
1640  if(edgeOrient==eBackwards)
1641  {
1642  reverseOrdering = true;
1643  }
1644  }
1645  else
1646  {
1647  IdxRange[0][0] = 2;
1648  IdxRange[0][1] = nummodes[0];
1649 
1650  if(edgeOrient==eBackwards)
1651  {
1652  signChange = true;
1653  }
1654  }
1655  }
1656  break;
1657  case 2:
1658  case 10:
1659  {
1660  if( bType[0] == LibUtilities::eGLL_Lagrange)
1661  {
1662  IdxRange[0][0] = 1;
1663  IdxRange[0][1] = nummodes[0] - 1;
1664 
1665  if(edgeOrient==eForwards)
1666  {
1667  reverseOrdering = true;
1668  }
1669  }
1670  else
1671  {
1672  IdxRange[0][0] = 2;
1673  IdxRange[0][1] = nummodes[0];
1674 
1675  if(edgeOrient==eForwards)
1676  {
1677  signChange = true;
1678  }
1679  }
1680  }
1681  break;
1682  }
1683 
1684  int p,q,r;
1685  int cnt = 0;
1686 
1687  for(r = IdxRange[2][0]; r < IdxRange[2][1]; r++)
1688  {
1689  for(q = IdxRange[1][0]; q < IdxRange[1][1]; q++)
1690  {
1691  for(p = IdxRange[0][0]; p < IdxRange[0][1]; p++)
1692  {
1693  maparray[cnt++]
1694  = r*nummodes[0]*nummodes[1] + q*nummodes[0] + p;
1695  }
1696  }
1697  }
1698 
1699  if( reverseOrdering )
1700  {
1701  reverse( maparray.get() , maparray.get()+nEdgeIntCoeffs );
1702  }
1703 
1704  if( signChange )
1705  {
1706  for(p = 1; p < nEdgeIntCoeffs; p+=2)
1707  {
1708  signarray[p] = -1;
1709  }
1710  }
1711  }
Principle Modified Functions .
Definition: BasisType.h:49
int GetEdgeNcoeffs(const int i) const
This function returns the number of expansion coefficients belonging to the i-th edge.
Definition: StdExpansion.h:287
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
Lagrange for SEM basis .
Definition: BasisType.h:53
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
Array< OneD, LibUtilities::BasisSharedPtr > m_base
int Nektar::StdRegions::StdHexExp::v_GetEdgeNcoeffs ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 740 of file StdHexExp.cpp.

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

741  {
742  ASSERTL2((i >= 0)&&(i <= 11),"edge id is out of range");
743 
744  if((i == 0)||(i == 2)||(i == 8)||(i == 10))
745  {
746  return GetBasisNumModes(0);
747  }
748  else if((i == 1)||(i == 3)||(i == 9)||(i == 11))
749  {
750  return GetBasisNumModes(1);
751  }
752  else
753  {
754  return GetBasisNumModes(2);
755  }
756  }
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:178
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:240
void Nektar::StdRegions::StdHexExp::v_GetFaceInteriorMap ( const int  fid,
const Orientation  faceOrient,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray 
)
protectedvirtual

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1718 of file StdHexExp.cpp.

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

1722  {
1725  "BasisType is not a boundary interior form");
1728  "BasisType is not a boundary interior form");
1731  "BasisType is not a boundary interior form");
1732 
1733  ASSERTL1((fid>=0)&&(fid<6),
1734  "local face id must be between 0 and 5");
1735 
1736  int nFaceIntCoeffs = GetFaceIntNcoeffs(fid);
1737 
1738  if(maparray.num_elements()!=nFaceIntCoeffs)
1739  {
1740  maparray = Array<OneD, unsigned int>(nFaceIntCoeffs);
1741  }
1742 
1743  if(signarray.num_elements() != nFaceIntCoeffs)
1744  {
1745  signarray = Array<OneD, int>(nFaceIntCoeffs,1);
1746  }
1747  else
1748  {
1749  fill( signarray.get() , signarray.get()+nFaceIntCoeffs, 1 );
1750  }
1751 
1752  int nummodes [3] = {m_base[0]->GetNumModes(),
1753  m_base[1]->GetNumModes(),
1754  m_base[2]->GetNumModes()};
1755 
1756  const LibUtilities::BasisType bType [3] = {GetBasisType(0),
1757  GetBasisType(1),
1758  GetBasisType(2)};
1759 
1760  int nummodesA = 0;
1761  int nummodesB = 0;
1762 
1763  // Determine the number of modes in face directions A & B based
1764  // on the face index given.
1765  switch(fid)
1766  {
1767  case 0:
1768  case 5:
1769  {
1770  nummodesA = nummodes[0];
1771  nummodesB = nummodes[1];
1772  }
1773  break;
1774  case 1:
1775  case 3:
1776  {
1777  nummodesA = nummodes[0];
1778  nummodesB = nummodes[2];
1779  }
1780  break;
1781  case 2:
1782  case 4:
1783  {
1784  nummodesA = nummodes[1];
1785  nummodesB = nummodes[2];
1786  }
1787  }
1788 
1789  int i,j;
1790  Array<OneD, int> arrayindx(nFaceIntCoeffs);
1791 
1792  // Create a mapping array to account for transposition of the
1793  // coordinates due to face orientation.
1794  for(i = 0; i < (nummodesB-2); i++)
1795  {
1796  for(j = 0; j < (nummodesA-2); j++)
1797  {
1798  if( faceOrient < 9 )
1799  {
1800  arrayindx[i*(nummodesA-2)+j] = i*(nummodesA-2)+j;
1801  }
1802  else
1803  {
1804  arrayindx[i*(nummodesA-2)+j] = j*(nummodesB-2)+i;
1805  }
1806  }
1807  }
1808 
1809  int IdxRange [3][2];
1810  int Incr[3];
1811 
1812  Array<OneD, int> sign0(nummodes[0], 1);
1813  Array<OneD, int> sign1(nummodes[1], 1);
1814  Array<OneD, int> sign2(nummodes[2], 1);
1815 
1816  // Set the upper and lower bounds, and increment for the faces
1817  // involving the first coordinate direction.
1818  switch(fid)
1819  {
1820  case 0: // bottom face
1821  {
1822  IdxRange[2][0] = 0;
1823  IdxRange[2][1] = 1;
1824  Incr[2] = 1;
1825  }
1826  break;
1827  case 5: // top face
1828  {
1829  if( bType[2] == LibUtilities::eGLL_Lagrange)
1830  {
1831  IdxRange[2][0] = nummodes[2] - 1;
1832  IdxRange[2][1] = nummodes[2];
1833  Incr[2] = 1;
1834  }
1835  else
1836  {
1837  IdxRange[2][0] = 1;
1838  IdxRange[2][1] = 2;
1839  Incr[2] = 1;
1840  }
1841 
1842  }
1843  break;
1844  default: // all other faces
1845  {
1846  if( bType[2] == LibUtilities::eGLL_Lagrange)
1847  {
1848  if( (((int) faceOrient)-5) % 2 )
1849  {
1850  IdxRange[2][0] = nummodes[2] - 2;
1851  IdxRange[2][1] = 0;
1852  Incr[2] = -1;
1853 
1854  }
1855  else
1856  {
1857  IdxRange[2][0] = 1;
1858  IdxRange[2][1] = nummodes[2] - 1;
1859  Incr[2] = 1;
1860  }
1861  }
1862  else
1863  {
1864  IdxRange[2][0] = 2;
1865  IdxRange[2][1] = nummodes[2];
1866  Incr[2] = 1;
1867 
1868  if( (((int) faceOrient)-5) % 2 )
1869  {
1870  for(i = 3; i < nummodes[2]; i+=2)
1871  {
1872  sign2[i] = -1;
1873  }
1874  }
1875  }
1876  }
1877  }
1878 
1879  // Set the upper and lower bounds, and increment for the faces
1880  // involving the second coordinate direction.
1881  switch(fid)
1882  {
1883  case 1:
1884  {
1885  IdxRange[1][0] = 0;
1886  IdxRange[1][1] = 1;
1887  Incr[1] = 1;
1888  }
1889  break;
1890  case 3:
1891  {
1892  if( bType[1] == LibUtilities::eGLL_Lagrange)
1893  {
1894  IdxRange[1][0] = nummodes[1] - 1;
1895  IdxRange[1][1] = nummodes[1];
1896  Incr[1] = 1;
1897  }
1898  else
1899  {
1900  IdxRange[1][0] = 1;
1901  IdxRange[1][1] = 2;
1902  Incr[1] = 1;
1903  }
1904  }
1905  break;
1906  case 0:
1907  case 5:
1908  {
1909  if( bType[1] == LibUtilities::eGLL_Lagrange)
1910  {
1911  if( (((int) faceOrient)-5) % 2 )
1912  {
1913  IdxRange[1][0] = nummodes[1] - 2;
1914  IdxRange[1][1] = 0;
1915  Incr[1] = -1;
1916 
1917  }
1918  else
1919  {
1920  IdxRange[1][0] = 1;
1921  IdxRange[1][1] = nummodes[1] - 1;
1922  Incr[1] = 1;
1923  }
1924  }
1925  else
1926  {
1927  IdxRange[1][0] = 2;
1928  IdxRange[1][1] = nummodes[1];
1929  Incr[1] = 1;
1930 
1931  if( (((int) faceOrient)-5) % 2 )
1932  {
1933  for(i = 3; i < nummodes[1]; i+=2)
1934  {
1935  sign1[i] = -1;
1936  }
1937  }
1938  }
1939  }
1940  break;
1941  default: // case2: case4:
1942  {
1943  if( bType[1] == LibUtilities::eGLL_Lagrange)
1944  {
1945  if( (((int) faceOrient)-5) % 4 > 1 )
1946  {
1947  IdxRange[1][0] = nummodes[1] - 2;
1948  IdxRange[1][1] = 0;
1949  Incr[1] = -1;
1950 
1951  }
1952  else
1953  {
1954  IdxRange[1][0] = 1;
1955  IdxRange[1][1] = nummodes[1] - 1;
1956  Incr[1] = 1;
1957  }
1958  }
1959  else
1960  {
1961  IdxRange[1][0] = 2;
1962  IdxRange[1][1] = nummodes[1];
1963  Incr[1] = 1;
1964 
1965  if( (((int) faceOrient)-5) % 4 > 1 )
1966  {
1967  for(i = 3; i < nummodes[1]; i+=2)
1968  {
1969  sign1[i] = -1;
1970  }
1971  }
1972  }
1973  }
1974  }
1975 
1976  switch(fid)
1977  {
1978  case 4:
1979  {
1980  IdxRange[0][0] = 0;
1981  IdxRange[0][1] = 1;
1982  Incr[0] = 1;
1983  }
1984  break;
1985  case 2:
1986  {
1987  if( bType[0] == LibUtilities::eGLL_Lagrange)
1988  {
1989  IdxRange[0][0] = nummodes[0] - 1;
1990  IdxRange[0][1] = nummodes[0];
1991  Incr[0] = 1;
1992  }
1993  else
1994  {
1995  IdxRange[0][0] = 1;
1996  IdxRange[0][1] = 2;
1997  Incr[0] = 1;
1998  }
1999  }
2000  break;
2001  default:
2002  {
2003  if( bType[0] == LibUtilities::eGLL_Lagrange)
2004  {
2005  if( (((int) faceOrient)-5) % 4 > 1 )
2006  {
2007  IdxRange[0][0] = nummodes[0] - 2;
2008  IdxRange[0][1] = 0;
2009  Incr[0] = -1;
2010 
2011  }
2012  else
2013  {
2014  IdxRange[0][0] = 1;
2015  IdxRange[0][1] = nummodes[0] - 1;
2016  Incr[0] = 1;
2017  }
2018  }
2019  else
2020  {
2021  IdxRange[0][0] = 2;
2022  IdxRange[0][1] = nummodes[0];
2023  Incr[0] = 1;
2024 
2025  if( (((int) faceOrient)-5) % 4 > 1 )
2026  {
2027  for(i = 3; i < nummodes[0]; i+=2)
2028  {
2029  sign0[i] = -1;
2030  }
2031  }
2032  }
2033  }
2034  }
2035 
2036  int p,q,r;
2037  int cnt = 0;
2038 
2039  for(r = IdxRange[2][0]; r != IdxRange[2][1]; r+=Incr[2])
2040  {
2041  for(q = IdxRange[1][0]; q != IdxRange[1][1]; q+=Incr[1])
2042  {
2043  for(p = IdxRange[0][0]; p != IdxRange[0][1]; p+=Incr[0])
2044  {
2045  maparray [ arrayindx[cnt ] ]
2046  = r*nummodes[0]*nummodes[1] + q*nummodes[0] + p;
2047  signarray[ arrayindx[cnt++] ]
2048  = sign0[p] * sign1[q] * sign2[r];
2049  }
2050  }
2051  }
2052  }
Principle Modified Functions .
Definition: BasisType.h:49
int GetFaceIntNcoeffs(const int i) const
Definition: StdExpansion.h:359
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
Lagrange for SEM basis .
Definition: BasisType.h:53
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
Array< OneD, LibUtilities::BasisSharedPtr > m_base
int Nektar::StdRegions::StdHexExp::v_GetFaceIntNcoeffs ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 782 of file StdHexExp.cpp.

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

783  {
784  ASSERTL2((i >= 0) && (i <= 5), "face id is out of range");
785  if((i == 0) || (i == 5))
786  {
787  return (GetBasisNumModes(0)-2)*(GetBasisNumModes(1)-2);
788  }
789  else if((i == 1) || (i == 3))
790  {
791  return (GetBasisNumModes(0)-2)*(GetBasisNumModes(2)-2);
792  }
793  else
794  {
795  return (GetBasisNumModes(1)-2)*(GetBasisNumModes(2)-2);
796  }
797 
798  }
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:178
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:240
int Nektar::StdRegions::StdHexExp::v_GetFaceNcoeffs ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 764 of file StdHexExp.cpp.

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

765  {
766  ASSERTL2((i >= 0) && (i <= 5), "face id is out of range");
767  if((i == 0) || (i == 5))
768  {
769  return GetBasisNumModes(0)*GetBasisNumModes(1);
770  }
771  else if((i == 1) || (i == 3))
772  {
773  return GetBasisNumModes(0)*GetBasisNumModes(2);
774  }
775  else
776  {
777  return GetBasisNumModes(1)*GetBasisNumModes(2);
778  }
779  }
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:178
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:240
int Nektar::StdRegions::StdHexExp::v_GetFaceNumPoints ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 807 of file StdHexExp.cpp.

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

808  {
809  ASSERTL2(i >= 0 && i <= 5, "face id is out of range");
810 
811  if (i == 0 || i == 5)
812  {
813  return m_base[0]->GetNumPoints()*
814  m_base[1]->GetNumPoints();
815  }
816  else if (i == 1 || i == 3)
817  {
818  return m_base[0]->GetNumPoints()*
819  m_base[2]->GetNumPoints();
820  }
821  else
822  {
823  return m_base[1]->GetNumPoints()*
824  m_base[2]->GetNumPoints();
825  }
826  }
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:240
Array< OneD, LibUtilities::BasisSharedPtr > m_base
LibUtilities::PointsKey Nektar::StdRegions::StdHexExp::v_GetFacePointsKey ( const int  i,
const int  j 
) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 828 of file StdHexExp.cpp.

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

830  {
831  ASSERTL2(i >= 0 && i <= 5, "face id is out of range");
832  ASSERTL2(j == 0 || j == 1, "face direction is out of range");
833 
834  if (i == 0 || i == 5)
835  {
836  return m_base[j]->GetPointsKey();
837  }
838  else if (i == 1 || i == 3)
839  {
840  return m_base[2*j]->GetPointsKey();
841  }
842  else
843  {
844  return m_base[j+1]->GetPointsKey();
845  }
846  }
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:240
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdHexExp::v_GetFaceToElementMap ( const int  fid,
const Orientation  faceOrient,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray,
int  P = -1,
int  Q = -1 
)
protectedvirtual

Only for basis type Modified_A or GLL_LAGRANGE in all directions.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 934 of file StdHexExp.cpp.

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

941  {
942  int i,j;
943  int nummodesA, nummodesB;
944 
947  "Method only implemented if BasisType is indentical in "
948  "all directions");
951  "Method only implemented for Modified_A or GLL_Lagrange BasisType");
952 
953  const int nummodes0 = m_base[0]->GetNumModes();
954  const int nummodes1 = m_base[1]->GetNumModes();
955  const int nummodes2 = m_base[2]->GetNumModes();
956 
957  switch(fid)
958  {
959  case 0:
960  case 5:
961  nummodesA = nummodes0;
962  nummodesB = nummodes1;
963  break;
964  case 1:
965  case 3:
966  nummodesA = nummodes0;
967  nummodesB = nummodes2;
968  break;
969  case 2:
970  case 4:
971  nummodesA = nummodes1;
972  nummodesB = nummodes2;
973  break;
974  }
975 
976  bool CheckForZeroedModes = false;
977 
978  if (P == -1)
979  {
980  P = nummodesA;
981  Q = nummodesB;
982  }
983 
984  if((P != nummodesA)||(Q != nummodesB))
985  {
986  CheckForZeroedModes = true;
987  }
988 
989  bool modified = (GetEdgeBasisType(0) == LibUtilities::eModified_A);
990  int nFaceCoeffs = P*Q;
991 
992  if(maparray.num_elements() != nFaceCoeffs)
993  {
994  maparray = Array<OneD, unsigned int>(nFaceCoeffs);
995  }
996 
997  if(signarray.num_elements() != nFaceCoeffs)
998  {
999  signarray = Array<OneD, int>(nFaceCoeffs,1);
1000  }
1001  else
1002  {
1003  fill( signarray.get() , signarray.get()+nFaceCoeffs, 1 );
1004  }
1005 
1006  Array<OneD, int> arrayindx(nFaceCoeffs);
1007 
1008  for(i = 0; i < Q; i++)
1009  {
1010  for(j = 0; j < P; j++)
1011  {
1012  if( faceOrient < 9 )
1013  {
1014  arrayindx[i*P+j] = i*P+j;
1015  }
1016  else
1017  {
1018  arrayindx[i*P+j] = j*Q+i;
1019  }
1020  }
1021  }
1022 
1023  int offset = 0;
1024  int jump1 = 1;
1025  int jump2 = 1;
1026 
1027  switch(fid)
1028  {
1029  case 5:
1030  {
1031  if (modified)
1032  {
1033  offset = nummodes0*nummodes1;
1034  }
1035  else
1036  {
1037  offset = (nummodes2-1)*nummodes0*nummodes1;
1038  jump1 = nummodes0;
1039  }
1040  }
1041  case 0:
1042  {
1043  jump1 = nummodes0;
1044  }
1045  break;
1046  case 3:
1047  {
1048  if (modified)
1049  {
1050  offset = nummodes0;
1051  }
1052  else
1053  {
1054  offset = nummodes0*(nummodes1-1);
1055  jump1 = nummodes0*nummodes1;
1056  }
1057  }
1058  case 1:
1059  {
1060  jump1 = nummodes0*nummodes1;
1061  }
1062  break;
1063  case 2:
1064  {
1065  if (modified)
1066  {
1067  offset = 1;
1068  }
1069  else
1070  {
1071  offset = nummodes0-1;
1072  jump1 = nummodes0*nummodes1;
1073  jump2 = nummodes0;
1074 
1075  }
1076  }
1077  case 4:
1078  {
1079  jump1 = nummodes0*nummodes1;
1080  jump2 = nummodes0;
1081  }
1082  break;
1083  default:
1084  ASSERTL0(false,"fid must be between 0 and 5");
1085  }
1086 
1087  for(i = 0; i < Q; i++)
1088  {
1089  for(j = 0; j < P; j++)
1090  {
1091  maparray[ arrayindx[i*P+j] ]
1092  = i*jump1 + j*jump2 + offset;
1093  }
1094  }
1095 
1096 
1097  if(CheckForZeroedModes)
1098  {
1099  if(modified)
1100  {
1101  // zero signmap and set maparray to zero if elemental
1102  // modes are not as large as face modesl
1103  for(i = 0; i < nummodesB; i++)
1104  {
1105  for(j = nummodesA; j < P; j++)
1106  {
1107  signarray[arrayindx[i*P+j]] = 0.0;
1108  maparray[arrayindx[i*P+j]] = maparray[0];
1109  }
1110  }
1111 
1112  for(i = nummodesB; i < Q; i++)
1113  {
1114  for(j = 0; j < P; j++)
1115  {
1116  signarray[arrayindx[i*P+j]] = 0.0;
1117  maparray[arrayindx[i*P+j]] = maparray[0];
1118  }
1119  }
1120  }
1121  else
1122  {
1123  ASSERTL0(false,"Different trace space face dimention and element face dimention not possible for GLL-Lagrange bases");
1124  }
1125  }
1126 
1127  if( (faceOrient==6) || (faceOrient==8) ||
1128  (faceOrient==11) || (faceOrient==12) )
1129  {
1130  if(faceOrient<9)
1131  {
1132  if (modified)
1133  {
1134  for(i = 3; i < Q; i+=2)
1135  {
1136  for(j = 0; j < P; j++)
1137  {
1138  signarray[ arrayindx[i*P+j] ] *= -1;
1139  }
1140  }
1141 
1142  for(i = 0; i < P; i++)
1143  {
1144  swap( maparray[i] , maparray[i+P] );
1145  swap( signarray[i] , signarray[i+P] );
1146  }
1147 
1148  }
1149  else
1150  {
1151  for(i = 0; i < P; i++)
1152  {
1153  for(j = 0; j < Q/2; j++)
1154  {
1155  swap( maparray[i + j*P],
1156  maparray[i+P*Q
1157  -P -j*P] );
1158  swap( signarray[i + j*P],
1159  signarray[i+P*Q
1160  -P -j*P]);
1161  }
1162  }
1163  }
1164  }
1165  else
1166  {
1167  if (modified)
1168  {
1169  for(i = 0; i < Q; i++)
1170  {
1171  for(j = 3; j < P; j+=2)
1172  {
1173  signarray[ arrayindx[i*P+j] ] *= -1;
1174  }
1175  }
1176 
1177  for(i = 0; i < Q; i++)
1178  {
1179  swap( maparray[i] , maparray[i+Q] );
1180  swap( signarray[i] , signarray[i+Q] );
1181  }
1182 
1183  }
1184  else
1185  {
1186  for(i = 0; i < P; i++)
1187  {
1188  for(j = 0; j < Q/2; j++)
1189  {
1190  swap( maparray[i*Q + j],
1191  maparray[i*Q + Q -1 -j]);
1192  swap( signarray[i*Q + j],
1193  signarray[i*Q + Q -1 -j]);
1194  }
1195  }
1196  }
1197  }
1198  }
1199 
1200  if( (faceOrient==7) || (faceOrient==8) ||
1201  (faceOrient==10) || (faceOrient==12) )
1202  {
1203  if(faceOrient<9)
1204  {
1205  if (modified)
1206  {
1207  for(i = 0; i < Q; i++)
1208  {
1209  for(j = 3; j < P; j+=2)
1210  {
1211  signarray[ arrayindx[i*P+j] ] *= -1;
1212  }
1213  }
1214 
1215  for(i = 0; i < Q; i++)
1216  {
1217  swap( maparray[i*P],
1218  maparray[i*P+1]);
1219  swap( signarray[i*P],
1220  signarray[i*P+1]);
1221  }
1222  }
1223  else
1224  {
1225  for(i = 0; i < Q; i++)
1226  {
1227  for(j = 0; j < P/2; j++)
1228  {
1229  swap( maparray[i*P + j],
1230  maparray[i*P + P -1 -j]);
1231  swap( signarray[i*P + j],
1232  signarray[i*P + P -1 -j]);
1233  }
1234  }
1235  }
1236 
1237 
1238 
1239  }
1240  else
1241  {
1242  if (modified)
1243  {
1244  for(i = 3; i < Q; i+=2)
1245  {
1246  for(j = 0; j < P; j++)
1247  {
1248  signarray[ arrayindx[i*P+j] ] *= -1;
1249  }
1250  }
1251 
1252  for(i = 0; i < P; i++)
1253  {
1254  swap( maparray[i*Q],
1255  maparray[i*Q+1]);
1256  swap( signarray[i*Q],
1257  signarray[i*Q+1]);
1258  }
1259  }
1260  else
1261  {
1262  for(i = 0; i < Q; i++)
1263  {
1264  for(j = 0; j < P/2; j++)
1265  {
1266  swap( maparray[i + j*Q] ,
1267  maparray[i+P*Q - Q -j*Q] );
1268  swap( signarray[i + j*Q] ,
1269  signarray[i+P*Q - Q -j*Q] );
1270  }
1271  }
1272  }
1273  }
1274  }
1275  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
Principle Modified Functions .
Definition: BasisType.h:49
LibUtilities::BasisType GetEdgeBasisType(const int i) const
This function returns the type of expansion basis on the i-th edge.
Definition: StdExpansion.h:413
Lagrange for SEM basis .
Definition: BasisType.h:53
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdHexExp::v_GetInteriorMap ( Array< OneD, unsigned int > &  outarray)
protectedvirtual
Parameters
outarrayStorage area for computed map.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2058 of file StdHexExp.cpp.

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

2059  {
2062  "BasisType is not a boundary interior form");
2065  "BasisType is not a boundary interior form");
2068  "BasisType is not a boundary interior form");
2069 
2070  int i;
2071  int nummodes [3] = {m_base[0]->GetNumModes(),
2072  m_base[1]->GetNumModes(),
2073  m_base[2]->GetNumModes()};
2074 
2075  int nIntCoeffs = m_ncoeffs - NumBndryCoeffs();
2076 
2077  if(outarray.num_elements() != nIntCoeffs)
2078  {
2079  outarray = Array<OneD, unsigned int>(nIntCoeffs);
2080  }
2081 
2082  const LibUtilities::BasisType Btype [3] = {GetBasisType(0),
2083  GetBasisType(1),
2084  GetBasisType(2)};
2085 
2086  int p,q,r;
2087  int cnt = 0;
2088 
2089  int IntIdx [3][2];
2090 
2091  for(i = 0; i < 3; i++)
2092  {
2093  if( Btype[i] == LibUtilities::eModified_A)
2094  {
2095  IntIdx[i][0] = 2;
2096  IntIdx[i][1] = nummodes[i];
2097  }
2098  else
2099  {
2100  IntIdx[i][0] = 1;
2101  IntIdx[i][1] = nummodes[i]-1;
2102  }
2103  }
2104 
2105  for(r = IntIdx[2][0]; r < IntIdx[2][1]; r++)
2106  {
2107  for( q = IntIdx[1][0]; q < IntIdx[1][1]; q++)
2108  {
2109  for( p = IntIdx[0][0]; p < IntIdx[0][1]; p++)
2110  {
2111  outarray[cnt++] = r*nummodes[0]*nummodes[1] +
2112  q*nummodes[0] + p;
2113  }
2114  }
2115  }
2116  }
Principle Modified Functions .
Definition: BasisType.h:49
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
Lagrange for SEM basis .
Definition: BasisType.h:53
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
Array< OneD, LibUtilities::BasisSharedPtr > m_base
int Nektar::StdRegions::StdHexExp::v_GetNedges ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 681 of file StdHexExp.cpp.

682  {
683  return 12;
684  }
int Nektar::StdRegions::StdHexExp::v_GetNfaces ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 687 of file StdHexExp.cpp.

688  {
689  return 6;
690  }
int Nektar::StdRegions::StdHexExp::v_GetNverts ( ) const
protectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 675 of file StdHexExp.cpp.

676  {
677  return 8;
678  }
int Nektar::StdRegions::StdHexExp::v_GetTotalEdgeIntNcoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 758 of file StdHexExp.cpp.

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

759  {
761  }
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:178
int Nektar::StdRegions::StdHexExp::v_GetTotalFaceIntNcoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 800 of file StdHexExp.cpp.

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

801  {
802  return 2*((GetBasisNumModes(0)-2)*(GetBasisNumModes(1)-2)+
803  (GetBasisNumModes(0)-2)*(GetBasisNumModes(2)-2)+
804  (GetBasisNumModes(1)-2)*(GetBasisNumModes(2)-2));
805  }
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:178
int Nektar::StdRegions::StdHexExp::v_GetVertexMap ( int  localVertexId,
bool  useCoeffPacking = false 
)
protectedvirtual

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1286 of file StdHexExp.cpp.

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

1287  {
1290  "BasisType is not a boundary interior form");
1293  "BasisType is not a boundary interior form");
1296  "BasisType is not a boundary interior form");
1297 
1298  ASSERTL1((localVertexId>=0)&&(localVertexId<8),
1299  "local vertex id must be between 0 and 7");
1300 
1301  int p = 0;
1302  int q = 0;
1303  int r = 0;
1304 
1305  // Retrieve the number of modes in each dimension.
1306  int nummodes [3] = {m_base[0]->GetNumModes(),
1307  m_base[1]->GetNumModes(),
1308  m_base[2]->GetNumModes()};
1309 
1310  if(useCoeffPacking == true) // follow packing of coefficients i.e q,r,p
1311  {
1312  if(localVertexId > 3)
1313  {
1315  {
1316  r = nummodes[2]-1;
1317  }
1318  else
1319  {
1320  r = 1;
1321  }
1322  }
1323 
1324  switch(localVertexId % 4)
1325  {
1326  case 0:
1327  break;
1328  case 1:
1329  {
1331  {
1332  p = nummodes[0]-1;
1333  }
1334  else
1335  {
1336  p = 1;
1337  }
1338  }
1339  break;
1340  case 2:
1341  {
1343  {
1344  q = nummodes[1]-1;
1345  }
1346  else
1347  {
1348  q = 1;
1349  }
1350  }
1351  break;
1352  case 3:
1353  {
1355  {
1356  p = nummodes[0]-1;
1357  q = nummodes[1]-1;
1358  }
1359  else
1360  {
1361  p = 1;
1362  q = 1;
1363  }
1364  }
1365  break;
1366  }
1367  }
1368  else
1369  {
1370  // Right face (vertices 1,2,5,6)
1371  if( (localVertexId % 4) % 3 > 0 )
1372  {
1374  {
1375  p = nummodes[0]-1;
1376  }
1377  else
1378  {
1379  p = 1;
1380  }
1381  }
1382 
1383  // Back face (vertices 2,3,6,7)
1384  if( localVertexId % 4 > 1 )
1385  {
1387  {
1388  q = nummodes[1]-1;
1389  }
1390  else
1391  {
1392  q = 1;
1393  }
1394  }
1395 
1396  // Top face (vertices 4,5,6,7)
1397  if( localVertexId > 3)
1398  {
1400  {
1401  r = nummodes[2]-1;
1402  }
1403  else
1404  {
1405  r = 1;
1406  }
1407  }
1408  }
1409  // Compute the local number.
1410  return r*nummodes[0]*nummodes[1] + q*nummodes[0] + p;
1411  }
Principle Modified Functions .
Definition: BasisType.h:49
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
Lagrange for SEM basis .
Definition: BasisType.h:53
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdHexExp::v_HelmholtzMatrixOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 2262 of file StdHexExp.cpp.

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

2266  {
2267  StdHexExp::v_HelmholtzMatrixOp_MatFree(inarray,outarray,mkey);
2268  }
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
void Nektar::StdRegions::StdHexExp::v_IProductWRTBase ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

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

Parameters
inarray?
outarray?

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 366 of file StdHexExp.cpp.

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

369  {
370  if(m_base[0]->Collocation() &&
371  m_base[1]->Collocation() &&
372  m_base[2]->Collocation())
373  {
374  MultiplyByQuadratureMetric(inarray,outarray);
375  }
376  else
377  {
378  StdHexExp::v_IProductWRTBase_SumFac(inarray,outarray);
379  }
380  }
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:942
Array< OneD, LibUtilities::BasisSharedPtr > m_base
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multbyweights=true)
Definition: StdHexExp.cpp:399
void Nektar::StdRegions::StdHexExp::v_IProductWRTBase_MatOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Implementation of the local matrix inner product operation.

Definition at line 385 of file StdHexExp.cpp.

References Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::StdRegions::eIProductWRTBase, Nektar::StdRegions::StdExpansion::GetStdMatrix(), Nektar::StdRegions::StdExpansion::GetTotPoints(), and Nektar::StdRegions::StdExpansion::m_ncoeffs.

387  {
388  int nq = GetTotPoints();
389  StdMatrixKey iprodmatkey(eIProductWRTBase,DetShapeType(),*this);
390  DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
391 
392  Blas::Dgemv('N',m_ncoeffs,nq,1.0,iprodmat->GetPtr().get(),
393  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
394  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:470
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:700
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:141
void Nektar::StdRegions::StdHexExp::v_IProductWRTBase_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
bool  multiplybyweights = true 
)
protectedvirtual

Implementation of the sum-factorization inner product operation.

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 399 of file StdHexExp.cpp.

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

Referenced by v_IProductWRTBase().

403  {
404  int nquad0 = m_base[0]->GetNumPoints();
405  int nquad1 = m_base[1]->GetNumPoints();
406  int nquad2 = m_base[2]->GetNumPoints();
407  int order0 = m_base[0]->GetNumModes();
408  int order1 = m_base[1]->GetNumModes();
409 
410  Array<OneD, NekDouble> wsp(nquad0*nquad1*(nquad2+order0) +
411  order0*order1*nquad2);
412 
413  if(multiplybyweights)
414  {
415  Array<OneD, NekDouble> tmp(inarray.num_elements());
416  MultiplyByQuadratureMetric(inarray,tmp);
417 
419  m_base[1]->GetBdata(),
420  m_base[2]->GetBdata(),
421  tmp,outarray,wsp,true,true,true);
422  }
423  else
424  {
426  m_base[1]->GetBdata(),
427  m_base[2]->GetBdata(),
428  inarray,outarray,wsp,true,true,true);
429  }
430  }
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:942
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)
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdHexExp::v_IProductWRTBase_SumFacKernel ( const Array< OneD, const NekDouble > &  base0,
const Array< OneD, const NekDouble > &  base1,
const Array< OneD, const NekDouble > &  base2,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
Array< OneD, NekDouble > &  wsp,
bool  doCheckCollDir0,
bool  doCheckCollDir1,
bool  doCheckCollDir2 
)
protectedvirtual

Implementation of the sum-factorisation inner product operation.

Todo:
Implement cases where only some directions are collocated.

Implements Nektar::StdRegions::StdExpansion3D.

Definition at line 437 of file StdHexExp.cpp.

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

446  {
447  int nquad0 = m_base[0]->GetNumPoints();
448  int nquad1 = m_base[1]->GetNumPoints();
449  int nquad2 = m_base[2]->GetNumPoints();
450  int nmodes0 = m_base[0]->GetNumModes();
451  int nmodes1 = m_base[1]->GetNumModes();
452  int nmodes2 = m_base[2]->GetNumModes();
453 
454  bool colldir0 = doCheckCollDir0?(m_base[0]->Collocation()):false;
455  bool colldir1 = doCheckCollDir1?(m_base[1]->Collocation()):false;
456  bool colldir2 = doCheckCollDir2?(m_base[2]->Collocation()):false;
457 
458  if(colldir0 && colldir1 && colldir2)
459  {
460  Vmath::Vcopy(m_ncoeffs,inarray.get(),1,outarray.get(),1);
461  }
462  else
463  {
464  ASSERTL1(wsp.num_elements() >= nmodes0*nquad2*(nquad1+nmodes1),
465  "Insufficient workspace size");
466 
467  Array<OneD, NekDouble> tmp0 = wsp;
468  Array<OneD, NekDouble> tmp1 = wsp + nmodes0*nquad1*nquad2;
469 
470 
471  if(colldir0)
472  {
473  // reshuffle data for next operation.
474  for(int n = 0; n < nmodes0; ++n)
475  {
476  Vmath::Vcopy(nquad1*nquad2,inarray.get()+n,nquad0,
477  tmp0.get()+nquad1*nquad2*n,1);
478  }
479  }
480  else
481  {
482  Blas::Dgemm('T', 'N', nquad1*nquad2, nmodes0, nquad0,
483  1.0, inarray.get(), nquad0,
484  base0.get(), nquad0,
485  0.0, tmp0.get(), nquad1*nquad2);
486  }
487 
488  if(colldir1)
489  {
490  // reshuffle data for next operation.
491  for(int n = 0; n < nmodes1; ++n)
492  {
493  Vmath::Vcopy(nquad2*nmodes0,tmp0.get()+n,nquad1,
494  tmp1.get()+nquad2*nmodes0*n,1);
495  }
496  }
497  else
498  {
499  Blas::Dgemm('T', 'N', nquad2*nmodes0, nmodes1, nquad1,
500  1.0, tmp0.get(), nquad1,
501  base1.get(), nquad1,
502  0.0, tmp1.get(), nquad2*nmodes0);
503  }
504 
505  if(colldir2)
506  {
507  // reshuffle data for next operation.
508  for(int n = 0; n < nmodes2; ++n)
509  {
510  Vmath::Vcopy(nmodes0*nmodes1,tmp1.get()+n,nquad2,
511  outarray.get()+nmodes0*nmodes1*n,1);
512  }
513  }
514  else
515  {
516  Blas::Dgemm('T', 'N', nmodes0*nmodes1, nmodes2, nquad2,
517  1.0, tmp1.get(), nquad2,
518  base2.get(), nquad2,
519  0.0, outarray.get(), nmodes0*nmodes1);
520  }
521  }
522  }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
void Nektar::StdRegions::StdHexExp::v_IProductWRTDerivBase ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 524 of file StdHexExp.cpp.

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

527  {
528  StdHexExp::IProductWRTDerivBase_SumFac(dir,inarray,outarray);
529  }
void IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void Nektar::StdRegions::StdHexExp::v_IProductWRTDerivBase_MatOp ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Definition at line 532 of file StdHexExp.cpp.

References ASSERTL0, Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::StdRegions::eIProductWRTDerivBase0, Nektar::StdRegions::eIProductWRTDerivBase1, Nektar::StdRegions::eIProductWRTDerivBase2, Nektar::StdRegions::StdExpansion::GetStdMatrix(), Nektar::StdRegions::StdExpansion::GetTotPoints(), and Nektar::StdRegions::StdExpansion::m_ncoeffs.

535  {
536  ASSERTL0((dir==0)||(dir==1)||(dir==2),"input dir is out of range");
537 
538  int nq = GetTotPoints();
539  MatrixType mtype;
540 
541  switch (dir)
542  {
543  case 0:
544  mtype = eIProductWRTDerivBase0;
545  break;
546  case 1:
547  mtype = eIProductWRTDerivBase1;
548  break;
549  case 2:
550  mtype = eIProductWRTDerivBase2;
551  break;
552  }
553 
554  StdMatrixKey iprodmatkey(mtype,DetShapeType(),*this);
555  DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
556 
557  Blas::Dgemv('N',m_ncoeffs,nq,1.0,iprodmat->GetPtr().get(),
558  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
559  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:470
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:700
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:141
void Nektar::StdRegions::StdHexExp::v_IProductWRTDerivBase_SumFac ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 562 of file StdHexExp.cpp.

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

565  {
566  ASSERTL0((dir==0)||(dir==1)||(dir==2),"input dir is out of range");
567 
568  int nquad1 = m_base[1]->GetNumPoints();
569  int nquad2 = m_base[2]->GetNumPoints();
570  int order0 = m_base[0]->GetNumModes();
571  int order1 = m_base[1]->GetNumModes();
572 
573  // If outarray > inarray then no need for temporary storage.
574  Array<OneD, NekDouble> tmp = outarray;
575  if (outarray.num_elements() < inarray.num_elements())
576  {
577  tmp = Array<OneD, NekDouble>(inarray.num_elements());
578  }
579 
580  // Need workspace for sumfackernel though
581  Array<OneD, NekDouble> wsp(order0*nquad2*(nquad1+order1));
582 
583  // multiply by integration constants
584  MultiplyByQuadratureMetric(inarray,tmp);
585 
586  // perform sum-factorisation
587  switch (dir)
588  {
589  case 0:
590  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
591  m_base[1]->GetBdata(),
592  m_base[2]->GetBdata(),
593  tmp,outarray,wsp,
594  false,true,true);
595  break;
596  case 1:
597  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
598  m_base[1]->GetDbdata(),
599  m_base[2]->GetBdata(),
600  tmp,outarray,wsp,
601  true,false,true);
602  break;
603  case 2:
604  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
605  m_base[1]->GetBdata(),
606  m_base[2]->GetDbdata(),
607  tmp,outarray,wsp,
608  true,true,false);
609  break;
610  }
611  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:942
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)
Array< OneD, LibUtilities::BasisSharedPtr > m_base
bool Nektar::StdRegions::StdHexExp::v_IsBoundaryInteriorExpansion ( )
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 85 of file StdHexExp.cpp.

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

86  {
87  return (m_base[0]->GetBasisType() == LibUtilities::eModified_A) &&
88  (m_base[1]->GetBasisType() == LibUtilities::eModified_A) &&
90  }
Principle Modified Functions .
Definition: BasisType.h:49
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdHexExp::v_LaplacianMatrixOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 2234 of file StdHexExp.cpp.

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

2238  {
2239  StdHexExp::v_LaplacianMatrixOp_MatFree(inarray,outarray,mkey);
2240  }
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
void Nektar::StdRegions::StdHexExp::v_LaplacianMatrixOp ( const int  k1,
const int  k2,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 2243 of file StdHexExp.cpp.

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

2247  {
2248  StdExpansion::LaplacianMatrixOp_MatFree(k1,k2,inarray,outarray,
2249  mkey);
2250  }
void LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void Nektar::StdRegions::StdHexExp::v_LocCoordToLocCollapsed ( const Array< OneD, const NekDouble > &  xi,
Array< OneD, NekDouble > &  eta 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 613 of file StdHexExp.cpp.

615  {
616  eta[0] = xi[0];
617  eta[1] = xi[1];
618  eta[2] = xi[2];
619  }
void Nektar::StdRegions::StdHexExp::v_MassMatrixOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 2225 of file StdHexExp.cpp.

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

2229  {
2230  StdExpansion::MassMatrixOp_MatFree(inarray,outarray,mkey);
2231  }
void MassMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void Nektar::StdRegions::StdHexExp::v_MultiplyByStdQuadratureMetric ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2294 of file StdHexExp.cpp.

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

2296  {
2297  int i;
2298  int nquad0 = m_base[0]->GetNumPoints();
2299  int nquad1 = m_base[1]->GetNumPoints();
2300  int nquad2 = m_base[2]->GetNumPoints();
2301  int nq01 = nquad0*nquad1;
2302  int nq12 = nquad1*nquad2;
2303 
2304  const Array<OneD, const NekDouble>& w0 = m_base[0]->GetW();
2305  const Array<OneD, const NekDouble>& w1 = m_base[1]->GetW();
2306  const Array<OneD, const NekDouble>& w2 = m_base[2]->GetW();
2307 
2308  for(i = 0; i < nq12; ++i)
2309  {
2310  Vmath::Vmul(nquad0, inarray.get()+i*nquad0, 1,
2311  w0.get(), 1, outarray.get()+i*nquad0,1);
2312  }
2313 
2314  for(i = 0; i < nq12; ++i)
2315  {
2316  Vmath::Smul(nquad0, w1[i%nquad1], outarray.get()+i*nquad0, 1,
2317  outarray.get()+i*nquad0, 1);
2318  }
2319 
2320  for(i = 0; i < nquad2; ++i)
2321  {
2322  Vmath::Smul(nq01, w2[i], outarray.get()+i*nq01, 1,
2323  outarray.get()+i*nq01, 1);
2324  }
2325  }
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:169
int Nektar::StdRegions::StdHexExp::v_NumBndryCoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 699 of file StdHexExp.cpp.

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

700  {
703  "BasisType is not a boundary interior form");
706  "BasisType is not a boundary interior form");
709  "BasisType is not a boundary interior form");
710 
711  int nmodes0 = m_base[0]->GetNumModes();
712  int nmodes1 = m_base[1]->GetNumModes();
713  int nmodes2 = m_base[2]->GetNumModes();
714 
715  return ( 2*( nmodes0*nmodes1 + nmodes0*nmodes2
716  + nmodes1*nmodes2)
717  - 4*( nmodes0 + nmodes1 + nmodes2 ) + 8 );
718  }
Principle Modified Functions .
Definition: BasisType.h:49
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
Lagrange for SEM basis .
Definition: BasisType.h:53
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
Array< OneD, LibUtilities::BasisSharedPtr > m_base
int Nektar::StdRegions::StdHexExp::v_NumDGBndryCoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 720 of file StdHexExp.cpp.

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

721  {
724  "BasisType is not a boundary interior form");
727  "BasisType is not a boundary interior form");
730  "BasisType is not a boundary interior form");
731 
732  int nmodes0 = m_base[0]->GetNumModes();
733  int nmodes1 = m_base[1]->GetNumModes();
734  int nmodes2 = m_base[2]->GetNumModes();
735 
736  return 2*( nmodes0*nmodes1 + nmodes0*nmodes2
737  + nmodes1*nmodes2 );
738  }
Principle Modified Functions .
Definition: BasisType.h:49
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
Lagrange for SEM basis .
Definition: BasisType.h:53
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdHexExp::v_PhysDeriv ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d0,
Array< OneD, NekDouble > &  out_d1,
Array< OneD, NekDouble > &  out_d2 
)
protectedvirtual

Differentiation Methods.

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 100 of file StdHexExp.cpp.

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

Referenced by v_StdPhysDeriv().

104  {
105  StdExpansion3D::PhysTensorDeriv(inarray, out_d0, out_d1, out_d2);
106  }
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d1, Array< OneD, NekDouble > &outarray_d2, Array< OneD, NekDouble > &outarray_d3)
Calculate the 3D derivative in the local tensor/collapsed coordinate at the physical points...
void Nektar::StdRegions::StdHexExp::v_PhysDeriv ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual
Parameters
dirDirection in which to compute derivative. Valid values are 0, 1, 2.
inarrayInput array.
outarrayOutput array.

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 115 of file StdHexExp.cpp.

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

118  {
119  switch(dir)
120  {
121  case 0:
122  {
123  PhysDeriv(inarray, outarray, NullNekDouble1DArray,
125  }
126  break;
127  case 1:
128  {
129  PhysDeriv(inarray, NullNekDouble1DArray, outarray,
131  }
132  break;
133  case 2:
134  {
136  NullNekDouble1DArray, outarray);
137  }
138  break;
139  default:
140  {
141  ASSERTL1(false,"input dir is out of range");
142  }
143  break;
144  }
145  }
static Array< OneD, NekDouble > NullNekDouble1DArray
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)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
void Nektar::StdRegions::StdHexExp::v_StdPhysDeriv ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d0,
Array< OneD, NekDouble > &  out_d1,
Array< OneD, NekDouble > &  out_d2 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 147 of file StdHexExp.cpp.

References v_PhysDeriv().

152  {
153  StdHexExp::v_PhysDeriv(inarray, out_d0, out_d1, out_d2);
154  }
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
Differentiation Methods.
Definition: StdHexExp.cpp:100
void Nektar::StdRegions::StdHexExp::v_StdPhysDeriv ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 157 of file StdHexExp.cpp.

References v_PhysDeriv().

160  {
161  StdHexExp::v_PhysDeriv(dir, inarray, outarray);
162  }
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
Differentiation Methods.
Definition: StdHexExp.cpp:100
void Nektar::StdRegions::StdHexExp::v_SVVLaplacianFilter ( Array< OneD, NekDouble > &  array,
const StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 2327 of file StdHexExp.cpp.

References Nektar::StdRegions::StdExpansion::BwdTrans(), Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::eFactorSVVDiffCoeff, Nektar::LibUtilities::eOrtho_A, Nektar::StdRegions::StdExpansion::FwdTrans(), Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), Nektar::StdRegions::StdExpansion::GetPointsType(), and Nektar::StdRegions::StdExpansion::m_base.

2329  {
2330  // Generate an orthonogal expansion
2331  int qa = m_base[0]->GetNumPoints();
2332  int qb = m_base[1]->GetNumPoints();
2333  int qc = m_base[2]->GetNumPoints();
2334  int nmodes_a = m_base[0]->GetNumModes();
2335  int nmodes_b = m_base[1]->GetNumModes();
2336  int nmodes_c = m_base[2]->GetNumModes();
2337  // Declare orthogonal basis.
2338  LibUtilities::PointsKey pa(qa,m_base[0]->GetPointsType());
2339  LibUtilities::PointsKey pb(qb,m_base[1]->GetPointsType());
2340  LibUtilities::PointsKey pc(qc,m_base[2]->GetPointsType());
2341 
2342  LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A,nmodes_a,pa);
2343  LibUtilities::BasisKey Bb(LibUtilities::eOrtho_A,nmodes_b,pb);
2344  LibUtilities::BasisKey Bc(LibUtilities::eOrtho_A,nmodes_c,pc);
2345  StdHexExp OrthoExp(Ba,Bb,Bc);
2346 
2347  Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
2348  int i,j,k;
2349 
2350  int cutoff = (int) (mkey.GetConstFactor(eFactorSVVCutoffRatio)*min(nmodes_a,nmodes_b));
2351  NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDiffCoeff);
2352 
2353  // project onto modal space.
2354  OrthoExp.FwdTrans(array,orthocoeffs);
2355 
2356 
2357  // Filter just trilinear space
2358  int nmodes = max(nmodes_a,nmodes_b);
2359  nmodes = max(nmodes,nmodes_c);
2360 
2361  Array<OneD, NekDouble> fac(nmodes,1.0);
2362  for(j = cutoff; j < nmodes; ++j)
2363  {
2364  fac[j] = fabs((j-nmodes)/((NekDouble) (j-cutoff+1.0)));
2365  fac[j] *= fac[j]; //added this line to conform with equation
2366  }
2367 
2368  for(i = 0; i < nmodes_a; ++i)
2369  {
2370  for(j = 0; j < nmodes_b; ++j)
2371  {
2372  for(k = 0; k < nmodes_c; ++k)
2373  {
2374  if((i >= cutoff)||(j >= cutoff)||(k >= cutoff))
2375  {
2376  orthocoeffs[i*nmodes_a*nmodes_b + j*nmodes_c + k] *= (SvvDiffCoeff*exp( -(fac[i]+fac[j]+fac[k]) ));
2377  }
2378  else
2379  {
2380  orthocoeffs[i*nmodes_a*nmodes_b + j*nmodes_c + k] *= 0.0;
2381  }
2382  }
2383  }
2384  }
2385 
2386  // backward transform to physical space
2387  OrthoExp.BwdTrans(orthocoeffs,array);
2388  }
Principle Orthogonal Functions .
Definition: BasisType.h:46
double NekDouble
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:216
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdHexExp::v_WeakDerivMatrixOp ( const int  i,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::HexExp.

Definition at line 2253 of file StdHexExp.cpp.

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

2257  {
2258  StdExpansion::WeakDerivMatrixOp_MatFree(i,inarray,outarray,
2259  mkey);
2260  }
void WeakDerivMatrixOp_MatFree(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)