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

#include <StdTetExp.h>

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

Public Member Functions

 StdTetExp ()
 StdTetExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
 StdTetExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc, NekDouble *coeffs, NekDouble *phys)
 StdTetExp (const StdTetExp &T)
 ~StdTetExp ()
LibUtilities::ShapeType DetShapeType () const
 This function returns the shape of the expansion domain.
NekDouble PhysEvaluate3D (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 Single Point Evaluation.
- 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.
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.
 StdExpansion (const int numcoeffs, const int numbases, const LibUtilities::BasisKey &Ba=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bb=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bc=LibUtilities::NullBasisKey)
 Constructor.
 StdExpansion (const StdExpansion &T)
 Copy Constructor.
virtual ~StdExpansion ()
 Destructor.
int GetNumBases () const
 This function returns the number of 1D bases used in the expansion.
const Array< OneD, const
LibUtilities::BasisSharedPtr > & 
GetBase () const
 This function gets the shared point to basis.
const
LibUtilities::BasisSharedPtr
GetBasis (int dir) const
 This function gets the shared point to basis in the dir direction.
int GetNcoeffs (void) const
 This function returns the total number of coefficients used in the expansion.
int GetTotPoints () const
 This function returns the total number of quadrature points used in the element.
LibUtilities::BasisType GetBasisType (const int dir) const
 This function returns the type of basis used in the dir direction.
int GetBasisNumModes (const int dir) const
 This function returns the number of expansion modes in the dir direction.
int EvalBasisNumModesMax (void) const
 This function returns the maximum number of expansion modes over all local directions.
LibUtilities::PointsType GetPointsType (const int dir) const
 This function returns the type of quadrature points used in the dir direction.
int GetNumPoints (const int dir) const
 This function returns the number of quadrature points in the dir direction.
const Array< OneD, const
NekDouble > & 
GetPoints (const int dir) const
 This function returns a pointer to the array containing the quadrature points in dir direction.
int GetNverts () const
 This function returns the number of vertices of the expansion domain.
int GetNedges () const
 This function returns the number of edges of the expansion domain.
int GetEdgeNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th edge.
int GetTotalEdgeIntNcoeffs () const
int GetEdgeNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th edge.
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.
int GetFaceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th face.
int GetFaceIntNcoeffs (const int i) const
int GetTotalFaceIntNcoeffs () const
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.
int GetNfaces () const
 This function returns the number of faces of the expansion domain.
int GetShapeDimension () const
bool IsBoundaryInteriorExpansion ()
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space.
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space.
void 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.
void FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 This function fills the array outarray with the mode-th mode of the expansion.
void IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 this function calculates the inner product of a given function f with the different modes of the expansion
void IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check)
void IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
int GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id.
void SetElmtId (const int id)
 Set the element id of this expansion when used in a list by returning value of m_elmt_id.
void GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2=NullNekDouble1DArray, Array< OneD, NekDouble > &coords_3=NullNekDouble1DArray)
 this function returns the physical coordinates of the quadrature points of the expansion
void GetCoord (const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
 given the coordinates of a point of the element in the local collapsed coordinate system, this function calculates the physical coordinates of the point
DNekMatSharedPtr GetStdMatrix (const StdMatrixKey &mkey)
DNekBlkMatSharedPtr GetStdStaticCondMatrix (const StdMatrixKey &mkey)
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, const Array< OneD, const 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)
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
void DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
StdRegions::Orientation GetFaceOrient (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)
void AddEdgeNormBoundaryInt (const int edge, boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
void AddEdgeNormBoundaryInt (const int edge, boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
void AddEdgeNormBoundaryBiInt (const int edge, boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
void AddFaceNormBoundaryInt (const int face, boost::shared_ptr< StdExpansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
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)
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.
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).
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 MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
DNekMatSharedPtr CreateGeneralMatrix (const StdMatrixKey &mkey)
 this function generates the mass matrix $\mathbf{M}[i][j] = \int \phi_i(\mathbf{x}) \phi_j(\mathbf{x}) d\mathbf{x}$
void GeneralMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
void 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)
void DGDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, boost::shared_ptr< StdExpansion > > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &coeffs, Array< OneD, NekDouble > &outarray)
NekDouble PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain.
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.
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta.
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.
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.
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 DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
virtual void v_DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
virtual StdRegions::Orientation v_GetFaceOrient (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.
NekDouble L2 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete $ L_2$ error, $ | \epsilon |_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 dx \right]^{1/2} d\xi_1 $ where $ u_{exact}$ is given by the array sol.
NekDouble H1 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete $ H^1$ error, $ | \epsilon |^1_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 + \nabla(u - u_{exact})\cdot\nabla(u - u_{exact})\cdot dx \right]^{1/2} d\xi_1 $ where $ u_{exact}$ is given by the array sol.
const 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)
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)
template<class T >
boost::shared_ptr< T > as ()

Protected Member Functions

virtual void v_PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dx, Array< OneD, NekDouble > &out_dy, Array< OneD, NekDouble > &out_dz)
 Calculate the derivative of the physical points.
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)
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_GetCoords (Array< OneD, NekDouble > &coords_x, Array< OneD, NekDouble > &coords_y, Array< OneD, NekDouble > &coords_z)
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_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_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.
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.
virtual void v_NegateFaceNormal (const int face)
- 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.
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.
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 IProductWRTBase_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)
virtual void v_AddEdgeNormBoundaryInt (const int edge, boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
virtual void v_AddEdgeNormBoundaryInt (const int edge, boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
virtual void v_AddFaceNormBoundaryInt (const int face, boost::shared_ptr< StdExpansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)

Private Member Functions

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

Additional Inherited Members

- Protected Attributes inherited from Nektar::StdRegions::StdExpansion3D
std::map< int, NormalVectorm_faceNormals
std::map< int, bool > m_negatedNormals

Detailed Description

Definition at line 51 of file StdTetExp.h.

Constructor & Destructor Documentation

Nektar::StdRegions::StdTetExp::StdTetExp ( )

Definition at line 45 of file StdTetExp.cpp.

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

Definition at line 50 of file StdTetExp.cpp.

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

:
Ba.GetNumModes(),
Bb.GetNumModes(),
Bc.GetNumModes()),
3, Ba, Bb, Bc),
Ba.GetNumModes(),
Bb.GetNumModes(),
Bc.GetNumModes()),
Ba, Bb, Bc)
{
ASSERTL0(Ba.GetNumModes() <= Bb.GetNumModes(),
"order in 'a' direction is higher than order "
"in 'b' direction");
ASSERTL0(Ba.GetNumModes() <= Bc.GetNumModes(),
"order in 'a' direction is higher than order "
"in 'c' direction");
ASSERTL0(Bb.GetNumModes() <= Bc.GetNumModes(),
"order in 'b' direction is higher than order "
"in 'c' direction");
}
Nektar::StdRegions::StdTetExp::StdTetExp ( const LibUtilities::BasisKey Ba,
const LibUtilities::BasisKey Bb,
const LibUtilities::BasisKey Bc,
NekDouble coeffs,
NekDouble phys 
)
Nektar::StdRegions::StdTetExp::StdTetExp ( const StdTetExp T)

Definition at line 75 of file StdTetExp.cpp.

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

Definition at line 82 of file StdTetExp.cpp.

{
}

Member Function Documentation

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

This function returns the shape of the expansion domain.

This function is a wrapper around the virtual function v_DetShapeType()

The different shape types implemented in the code are defined in the ::ShapeType enumeration list. As a result, the function will return one of the types of this enumeration list.

Returns
returns the shape of the expansion domain

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 70 of file StdTetExp.h.

References Nektar::LibUtilities::eTetrahedron.

Referenced by Nektar::LocalRegions::TetExp::CreateMatrix(), Nektar::StdRegions::StdNodalTetExp::ModalToNodal(), Nektar::StdRegions::StdNodalTetExp::NodalToModal(), Nektar::StdRegions::StdNodalTetExp::NodalToModalTranspose(), v_DetShapeType(), Nektar::LocalRegions::TetExp::v_FwdTrans(), Nektar::StdRegions::StdNodalTetExp::v_FwdTrans(), v_FwdTrans(), v_IProductWRTBase_MatOp(), and v_IProductWRTDerivBase_MatOp().

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

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

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

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

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

Note that in this element, we must have that $ P \leq Q \leq R$.

Definition at line 1752 of file StdTetExp.cpp.

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

Referenced by v_GetBoundaryMap(), v_GetEdgeInteriorMap(), v_GetFaceInteriorMap(), v_GetFaceToElementMap(), v_GetInteriorMap(), and v_GetVertexMap().

{
const int Q = m_base[1]->GetNumModes();
const int R = m_base[2]->GetNumModes();
int i,j,q_hat,k_hat;
int cnt = 0;
// Traverse to q-r plane number I
for (i = 0; i < I; ++i)
{
// Size of triangle part
q_hat = min(Q,R-i);
// Size of rectangle part
k_hat = max(R-Q-i,0);
cnt += q_hat*(q_hat+1)/2 + k_hat*Q;
}
// Traverse to q column J
q_hat = R-I;
for (j = 0; j < J; ++j)
{
cnt += q_hat;
q_hat--;
}
// Traverse up stacks to K
cnt += K;
return cnt;
}
NekDouble Nektar::StdRegions::StdTetExp::PhysEvaluate3D ( const Array< OneD, const NekDouble > &  coords,
const Array< OneD, const NekDouble > &  physvals 
)

Single Point Evaluation.

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

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

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

$ f_{pq} (\xi_{3k}) = \sum_{r=0}^{Q_z} \hat u_{pqr} \psi_{pqr}^c (\xi_{3k}),\\ g_{p} (\xi_{2j}, \xi_{3k}) = \sum_{r=0}^{Q_y} \psi_{pq}^b (\xi_{2j}) f_{pq} (\xi_{3k}),\\ u(\xi_{1i}, \xi_{2j}, \xi_{3k}) = \sum_{p=0}^{Q_x} \psi_{p}^a (\xi_{1i}) g_{p} (\xi_{2j}, \xi_{3k}). $

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTetExp.

Definition at line 293 of file StdTetExp.cpp.

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

Referenced by v_FillMode().

{
"Basis[1] is not a general tensor type");
"Basis[2] is not a general tensor type");
if(m_base[0]->Collocation() && m_base[1]->Collocation()
&& m_base[2]->Collocation())
{
inarray, 1, outarray, 1);
}
else
{
StdTetExp::v_BwdTrans_SumFac(inarray,outarray);
}
}
void Nektar::StdRegions::StdTetExp::v_BwdTrans_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Sum-factorisation implementation of the BwdTrans operation.

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTetExp.

Definition at line 323 of file StdTetExp.cpp.

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

Referenced by v_BwdTrans().

{
int nquad1 = m_base[1]->GetNumPoints();
int nquad2 = m_base[2]->GetNumPoints();
int order0 = m_base[0]->GetNumModes();
int order1 = m_base[1]->GetNumModes();
Array<OneD, NekDouble> wsp(nquad2*order0*order1*(order1+1)/2+
nquad2*nquad1*order0);
BwdTrans_SumFacKernel(m_base[0]->GetBdata(),
m_base[1]->GetBdata(),
m_base[2]->GetBdata(),
inarray,outarray,wsp,
true,true,true);
}
void Nektar::StdRegions::StdTetExp::v_BwdTrans_SumFacKernel ( const Array< OneD, const NekDouble > &  base0,
const Array< OneD, const NekDouble > &  base1,
const Array< OneD, const NekDouble > &  base2,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
Array< OneD, NekDouble > &  wsp,
bool  doCheckCollDir0,
bool  doCheckCollDir1,
bool  doCheckCollDir2 
)
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 356 of file StdTetExp.cpp.

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

{
int nquad0 = m_base[0]->GetNumPoints();
int nquad1 = m_base[1]->GetNumPoints();
int nquad2 = m_base[2]->GetNumPoints();
int order0 = m_base[0]->GetNumModes();
int order1 = m_base[1]->GetNumModes();
int order2 = m_base[2]->GetNumModes();
Array<OneD, NekDouble > tmp = wsp;
Array<OneD, NekDouble > tmp1 = tmp + nquad2*order0*order1*(order1+1)/2;
//Array<OneD, NekDouble > tmp(nquad2*order0*(order1+1)/2);
//Array<OneD, NekDouble > tmp1(nquad2*nquad1*order0);
int i, j, mode, mode1, cnt;
// Perform summation over '2' direction
mode = mode1 = cnt = 0;
for(i = 0; i < order0; ++i)
{
for(j = 0; j < order1-i; ++j, ++cnt)
{
Blas::Dgemv('N', nquad2, order2-i-j,
1.0, base2.get()+mode*nquad2, nquad2,
inarray.get()+mode1, 1,
0.0, tmp.get()+cnt*nquad2, 1);
mode += order2-i-j;
mode1 += order2-i-j;
}
//increment mode in case order1!=order2
for(j = order1-i; j < order2-i; ++j)
{
mode += order2-i-j;
}
}
// fix for modified basis by adding split of top singular
// vertex mode - currently (1+c)/2 x (1-b)/2 x (1-a)/2
// component is evaluated
{
// top singular vertex - (1+c)/2 x (1+b)/2 x (1-a)/2 component
Blas::Daxpy(nquad2,inarray[1],base2.get()+nquad2,1,
&tmp[0]+nquad2,1);
// top singular vertex - (1+c)/2 x (1-b)/2 x (1+a)/2 component
Blas::Daxpy(nquad2,inarray[1],base2.get()+nquad2,1,
&tmp[0]+order1*nquad2,1);
}
// Perform summation over '1' direction
mode = 0;
for(i = 0; i < order0; ++i)
{
Blas::Dgemm('N', 'T', nquad1, nquad2, order1-i,
1.0, base1.get()+mode*nquad1, nquad1,
tmp.get()+mode*nquad2, nquad2,
0.0, tmp1.get()+i*nquad1*nquad2, nquad1);
mode += order1-i;
}
// fix for modified basis by adding additional split of
// top and base singular vertex modes as well as singular
// edge
{
// use tmp to sort out singular vertices and
// singular edge components with (1+b)/2 (1+a)/2 form
for(i = 0; i < nquad2; ++i)
{
Blas::Daxpy(nquad1,tmp[nquad2+i], base1.get()+nquad1,1,
&tmp1[nquad1*nquad2]+i*nquad1,1);
}
}
// Perform summation over '0' direction
Blas::Dgemm('N', 'T', nquad0, nquad1*nquad2, order0,
1.0, base0.get(), nquad0,
tmp1.get(), nquad1*nquad2,
0.0, outarray.get(), nquad0);
}
int Nektar::StdRegions::StdTetExp::v_CalcNumberOfCoefficients ( const std::vector< unsigned int > &  nummodes,
int &  modes_offset 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1094 of file StdTetExp.cpp.

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

{
nummodes[modes_offset],
nummodes[modes_offset+1],
nummodes[modes_offset+2]);
modes_offset += 3;
return nmodes;
}
DNekMatSharedPtr Nektar::StdRegions::StdTetExp::v_CreateStdMatrix ( const StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TetExp, and Nektar::StdRegions::StdNodalTetExp.

Definition at line 1721 of file StdTetExp.cpp.

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

const LibUtilities::BasisKey Nektar::StdRegions::StdTetExp::v_DetFaceBasisKey ( const int  i,
const int  k 
) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1107 of file StdTetExp.cpp.

References ASSERTL2, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::StdRegions::StdExpansion::GetBasis(), and Nektar::LibUtilities::NullBasisKey().

{
ASSERTL2(i >= 0 && i <= 4, "face id is out of range");
ASSERTL2(k == 0 || k == 1, "face direction out of range");
int nummodes = GetBasis(0)->GetNumModes();
//temporary solution, need to add conditions based on face id
//also need to add check of the points type
switch (k)
{
case 0:
{
const LibUtilities::PointsKey pkey(nummodes+1,LibUtilities::eGaussLobattoLegendre);
return LibUtilities::BasisKey(LibUtilities::eModified_A,nummodes,pkey);
}
break;
case 1:
{
//const LibUtilities::PointsKey pkey(nummodes,LibUtilities::eGaussRadauMAlpha1Beta0);
const LibUtilities::PointsKey pkey(nummodes+1,LibUtilities::eGaussLobattoLegendre);
return LibUtilities::BasisKey(LibUtilities::eModified_B,nummodes,pkey);
}
break;
}
// Should not get here.
}
LibUtilities::ShapeType Nektar::StdRegions::StdTetExp::v_DetShapeType ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TetExp.

Definition at line 918 of file StdTetExp.cpp.

References DetShapeType().

{
return DetShapeType();
}
void Nektar::StdRegions::StdTetExp::v_FillMode ( const int  mode,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTetExp.

Definition at line 889 of file StdTetExp.cpp.

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

{
Array<OneD, NekDouble> tmp(m_ncoeffs,0.0);
tmp[mode] = 1.0;
StdTetExp::v_BwdTrans(tmp, outarray);
}
void Nektar::StdRegions::StdTetExp::v_FwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual
Parameters
inarrayarray of physical quadrature points to be transformed.
outarrayupdated array of expansion coefficients.

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTetExp, and Nektar::LocalRegions::TetExp.

Definition at line 455 of file StdTetExp.cpp.

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

{
v_IProductWRTBase(inarray,outarray);
// get Mass matrix inverse
StdMatrixKey masskey(eInvMass,DetShapeType(),*this);
DNekMatSharedPtr matsys = GetStdMatrix(masskey);
// copy inarray in case inarray == outarray
DNekVec in (m_ncoeffs, outarray);
DNekVec out(m_ncoeffs, outarray, eWrapper);
out = (*matsys)*in;
}
DNekMatSharedPtr Nektar::StdRegions::StdTetExp::v_GenMatrix ( const StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TetExp, and Nektar::StdRegions::StdNodalTetExp.

Definition at line 1716 of file StdTetExp.cpp.

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

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

List of all boundary modes in the the expansion.

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTetExp.

Definition at line 1663 of file StdTetExp.cpp.

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

{
"BasisType is not a boundary interior form");
"BasisType is not a boundary interior form");
"BasisType is not a boundary interior form");
int P = m_base[0]->GetNumModes();
int Q = m_base[1]->GetNumModes();
int R = m_base[2]->GetNumModes();
int i,j,k;
int idx = 0;
for (i = 0; i < P; ++i)
{
// First two Q-R planes are entirely boundary modes
if (i < 2)
{
for (j = 0; j < Q-i; j++)
{
for (k = 0; k < R-i-j; ++k)
{
outarray[idx++] = GetMode(i,j,k);
}
}
}
// Remaining Q-R planes contain boundary modes on bottom and
// left edge.
else
{
for (k = 0; k < R-i; ++k)
{
outarray[idx++] = GetMode(i,0,k);
}
for (j = 1; j < Q-i; ++j)
{
outarray[idx++] = GetMode(i,j,0);
}
}
}
}
void Nektar::StdRegions::StdTetExp::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::TetExp.

Definition at line 1154 of file StdTetExp.cpp.

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

{
Array<OneD, const NekDouble> eta_x = m_base[0]->GetZ();
Array<OneD, const NekDouble> eta_y = m_base[1]->GetZ();
Array<OneD, const NekDouble> eta_z = m_base[2]->GetZ();
int Qx = GetNumPoints(0);
int Qy = GetNumPoints(1);
int Qz = GetNumPoints(2);
// Convert collapsed coordinates into cartesian coordinates: eta
// --> xi
for( int k = 0; k < Qz; ++k ) {
for( int j = 0; j < Qy; ++j ) {
for( int i = 0; i < Qx; ++i ) {
int s = i + Qx*(j + Qy*k);
xi_x[s] = (eta_x[i] + 1.0) * (1.0 - eta_y[j]) * (1.0 - eta_z[k]) / 4 - 1.0;
xi_y[s] = (eta_y[j] + 1.0) * (1.0 - eta_z[k]) / 2 - 1.0;
xi_z[s] = eta_z[k];
}
}
}
}
LibUtilities::BasisType Nektar::StdRegions::StdTetExp::v_GetEdgeBasisType ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1136 of file StdTetExp.cpp.

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

{
ASSERTL2(i >= 0 && i <= 5, "edge id is out of range");
if (i == 0)
{
return GetBasisType(0);
}
else if (i == 1 || i == 2)
{
return GetBasisType(1);
}
else
{
return GetBasisType(2);
}
}
void Nektar::StdRegions::StdTetExp::v_GetEdgeInteriorMap ( const int  eid,
const Orientation  edgeOrient,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray 
)
protectedvirtual

Maps interior modes of an edge to the elemental modes.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1413 of file StdTetExp.cpp.

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

{
int i;
const int P = m_base[0]->GetNumModes();
const int Q = m_base[1]->GetNumModes();
const int R = m_base[2]->GetNumModes();
const int nEdgeIntCoeffs = v_GetEdgeNcoeffs(eid)-2;
if(maparray.num_elements() != nEdgeIntCoeffs)
{
maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
}
else
{
fill( maparray.get(), maparray.get() + nEdgeIntCoeffs, 0);
}
if(signarray.num_elements() != nEdgeIntCoeffs)
{
signarray = Array<OneD, int>(nEdgeIntCoeffs,1);
}
else
{
fill( signarray.get() , signarray.get()+nEdgeIntCoeffs, 1 );
}
switch (eid)
{
case 0:
for (i = 0; i < P-2; ++i)
{
maparray[i] = GetMode(i+2, 0, 0);
}
if(edgeOrient==eBackwards)
{
for(i = 1; i < nEdgeIntCoeffs; i+=2)
{
signarray[i] = -1;
}
}
break;
case 1:
for (i = 0; i < Q-2; ++i)
{
maparray[i] = GetMode(1, i+1, 0);
}
if(edgeOrient==eBackwards)
{
for(i = 1; i < nEdgeIntCoeffs; i+=2)
{
signarray[i] = -1;
}
}
break;
case 2:
for (i = 0; i < Q-2; ++i)
{
maparray[i] = GetMode(0, i+2, 0);
}
if(edgeOrient==eBackwards)
{
for(i = 1; i < nEdgeIntCoeffs; i+=2)
{
signarray[i] = -1;
}
}
break;
case 3:
for (i = 0; i < R-2; ++i)
{
maparray[i] = GetMode(0, 0, i+2);
}
if(edgeOrient==eBackwards)
{
for(i = 1; i < nEdgeIntCoeffs; i+=2)
{
signarray[i] = -1;
}
}
break;
case 4:
for (i = 0; i < R - 2; ++i)
{
maparray[i] = GetMode(1, 0, i+1);
}
if(edgeOrient==eBackwards)
{
for(i = 1; i < nEdgeIntCoeffs; i+=2)
{
signarray[i] = -1;
}
}
break;
case 5:
for (i = 0; i < R-2; ++i)
{
maparray[i] = GetMode(0, 1, i+1);
}
if(edgeOrient==eBackwards)
{
for(i = 1; i < nEdgeIntCoeffs; i+=2)
{
signarray[i] = -1;
}
}
break;
default:
ASSERTL0(false, "Edge not defined.");
break;
}
}
int Nektar::StdRegions::StdTetExp::v_GetEdgeNcoeffs ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 965 of file StdTetExp.cpp.

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

Referenced by v_GetEdgeInteriorMap().

{
ASSERTL2((i >= 0) && (i <= 5), "edge id is out of range");
int P = m_base[0]->GetNumModes();
int Q = m_base[1]->GetNumModes();
int R = m_base[2]->GetNumModes();
if (i == 0)
{
return P;
}
else if (i == 1 || i == 2)
{
return Q;
}
else
{
return R;
}
}
void Nektar::StdRegions::StdTetExp::v_GetFaceInteriorMap ( const int  fid,
const Orientation  faceOrient,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1530 of file StdTetExp.cpp.

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

{
int i, j, idx, k;
const int P = m_base[0]->GetNumModes();
const int Q = m_base[1]->GetNumModes();
const int R = m_base[2]->GetNumModes();
const int nFaceIntCoeffs = v_GetFaceIntNcoeffs(fid);
if(maparray.num_elements() != nFaceIntCoeffs)
{
maparray = Array<OneD, unsigned int>(nFaceIntCoeffs);
}
if(signarray.num_elements() != nFaceIntCoeffs)
{
signarray = Array<OneD, int>(nFaceIntCoeffs,1);
}
else
{
fill( signarray.get() , signarray.get()+nFaceIntCoeffs, 1 );
}
switch (fid)
{
case 0:
idx = 0;
for (i = 2; i < P-1; ++i)
{
for (j = 1; j < Q-i; ++j)
{
if ((int)faceOrient == 7)
{
signarray[idx] = (i%2 ? -1 : 1);
}
maparray[idx++] = GetMode(i,j,0);
}
}
break;
case 1:
idx = 0;
for (i = 2; i < P; ++i)
{
for (k = 1; k < R-i; ++k)
{
if ((int)faceOrient == 7)
{
signarray[idx] = (i%2 ? -1: 1);
}
maparray[idx++] = GetMode(i,0,k);
}
}
break;
case 2:
idx = 0;
for (j = 1; j < Q-2; ++j)
{
for (k = 1; k < R-1-j; ++k)
{
if ((int)faceOrient == 7)
{
signarray[idx] = ((j+1)%2 ? -1: 1);
}
maparray[idx++] = GetMode(1,j,k);
}
}
break;
case 3:
idx = 0;
for (j = 2; j < Q-1; ++j)
{
for (k = 1; k < R-j; ++k)
{
if ((int)faceOrient == 7)
{
signarray[idx] = (j%2 ? -1: 1);
}
maparray[idx++] = GetMode(0,j,k);
}
}
break;
default:
ASSERTL0(false, "Face interior map not available.");
break;
}
}
int Nektar::StdRegions::StdTetExp::v_GetFaceIntNcoeffs ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1021 of file StdTetExp.cpp.

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

Referenced by v_GetFaceInteriorMap().

{
ASSERTL2((i >= 0) && (i <= 3), "face id is out of range");
int Pi = m_base[0]->GetNumModes() - 2;
int Qi = m_base[1]->GetNumModes() - 2;
int Ri = m_base[2]->GetNumModes() - 2;
if((i == 0))
{
return Pi * (2*Qi - Pi - 1) / 2;
}
else if((i == 1))
{
return Pi * (2*Ri - Pi - 1) / 2;
}
else
{
return Qi * (2*Ri - Qi - 1) / 2;
}
}
int Nektar::StdRegions::StdTetExp::v_GetFaceNcoeffs ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 995 of file StdTetExp.cpp.

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

{
ASSERTL2((i >= 0) && (i <= 3), "face id is out of range");
int nFaceCoeffs = 0;
int nummodesA, nummodesB, P, Q;
if (i == 0)
{
nummodesA = GetBasisNumModes(0);
nummodesB = GetBasisNumModes(1);
}
else if ((i == 1) || (i == 2))
{
nummodesA = GetBasisNumModes(0);
nummodesB = GetBasisNumModes(2);
}
else
{
nummodesA = GetBasisNumModes(1);
nummodesB = GetBasisNumModes(2);
}
P = nummodesA - 1;
Q = nummodesB - 1;
nFaceCoeffs = Q+1 + (P*(1 + 2*Q - P))/2;
return nFaceCoeffs;
}
int Nektar::StdRegions::StdTetExp::v_GetFaceNumPoints ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1053 of file StdTetExp.cpp.

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

{
ASSERTL2(i >= 0 && i <= 3, "face id is out of range");
if (i == 0)
{
return m_base[0]->GetNumPoints()*
m_base[1]->GetNumPoints();
}
else if (i == 1)
{
return m_base[0]->GetNumPoints()*
m_base[2]->GetNumPoints();
}
else
{
return m_base[1]->GetNumPoints()*
m_base[2]->GetNumPoints();
}
}
LibUtilities::PointsKey Nektar::StdRegions::StdTetExp::v_GetFacePointsKey ( const int  i,
const int  j 
) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1074 of file StdTetExp.cpp.

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

{
ASSERTL2(i >= 0 && i <= 3, "face id is out of range");
ASSERTL2(j == 0 || j == 1, "face direction is out of range");
if (i == 0)
{
return m_base[j]->GetPointsKey();
}
else if (i == 1)
{
return m_base[2*j]->GetPointsKey();
}
else
{
return m_base[j+1]->GetPointsKey();
}
}
void Nektar::StdRegions::StdTetExp::v_GetFaceToElementMap ( const int  fid,
const Orientation  faceOrient,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray,
int  nummodesA = -1,
int  nummodesB = -1 
)
protectedvirtual

Maps Expansion2D modes of a 2D face to the corresponding expansion modes.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1196 of file StdTetExp.cpp.

References ASSERTL0, ASSERTL1, GetMode(), Nektar::StdRegions::StdExpansion::m_base, and v_IsBoundaryInteriorExpansion().

{
int P, Q, i, j, k, idx;
"Method only implemented for Modified_A BasisType (x "
"direction), Modified_B BasisType (y direction), and "
"Modified_C BasisType(z direction)");
int nFaceCoeffs = 0;
if (nummodesA == -1)
{
switch(fid)
{
case 0:
nummodesA = m_base[0]->GetNumModes();
nummodesB = m_base[1]->GetNumModes();
break;
case 1:
nummodesA = m_base[0]->GetNumModes();
nummodesB = m_base[2]->GetNumModes();
break;
case 2:
case 3:
nummodesA = m_base[1]->GetNumModes();
nummodesB = m_base[2]->GetNumModes();
break;
}
}
P = nummodesA;
Q = nummodesB;
nFaceCoeffs = Q + ((P-1)*(1 + 2*(Q-1) - (P-1)))/2;
// Allocate the map array and sign array; set sign array to ones (+)
if(maparray.num_elements() != nFaceCoeffs)
{
maparray = Array<OneD, unsigned int>(nFaceCoeffs,1);
}
if(signarray.num_elements() != nFaceCoeffs)
{
signarray = Array<OneD, int>(nFaceCoeffs,1);
}
else
{
fill( signarray.get() , signarray.get()+nFaceCoeffs, 1 );
}
switch (fid)
{
case 0:
idx = 0;
for (i = 0; i < P; ++i)
{
for (j = 0; j < Q-i; ++j)
{
if ((int)faceOrient == 7 && i > 1)
{
signarray[idx] = (i%2 ? -1 : 1);
}
maparray[idx++] = GetMode(i,j,0);
}
}
break;
case 1:
idx = 0;
for (i = 0; i < P; ++i)
{
for (k = 0; k < Q-i; ++k)
{
if ((int)faceOrient == 7 && i > 1)
{
signarray[idx] = (i%2 ? -1: 1);
}
maparray[idx++] = GetMode(i,0,k);
}
}
break;
case 2:
idx = 0;
for (j = 0; j < P-1; ++j)
{
for (k = 0; k < Q-1-j; ++k)
{
if ((int)faceOrient == 7 && j > 1)
{
signarray[idx] = ((j+1)%2 ? -1: 1);
}
maparray[idx++] = GetMode(1,j,k);
// Incorporate modes from zeroth plane where needed.
if (j == 0 && k == 0) {
maparray[idx++] = GetMode(0,0,1);
}
if (j == 0 && k == Q-2) {
for (int r = 0; r < Q-1; ++r) {
maparray[idx++] = GetMode(0,1,r);
}
}
}
}
break;
case 3:
idx = 0;
for (j = 0; j < P; ++j)
{
for (k = 0; k < Q-j; ++k)
{
if ((int)faceOrient == 7 && j > 1)
{
signarray[idx] = (j%2 ? -1: 1);
}
maparray[idx++] = GetMode(0,j,k);
}
}
break;
default:
ASSERTL0(false, "Element map not available.");
}
if ((int)faceOrient == 7)
{
swap(maparray[0], maparray[Q]);
for (i = 1; i < Q-1; ++i)
{
swap(maparray[i+1], maparray[Q+i]);
}
}
}
void Nektar::StdRegions::StdTetExp::v_GetInteriorMap ( Array< OneD, unsigned int > &  outarray)
protectedvirtual

List of all interior modes in the expansion.

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTetExp.

Definition at line 1624 of file StdTetExp.cpp.

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

{
"BasisType is not a boundary interior form");
"BasisType is not a boundary interior form");
"BasisType is not a boundary interior form");
int P = m_base[0]->GetNumModes();
int Q = m_base[1]->GetNumModes();
int R = m_base[2]->GetNumModes();
int nIntCoeffs = m_ncoeffs - NumBndryCoeffs();
if(outarray.num_elements() != nIntCoeffs)
{
outarray = Array<OneD, unsigned int>(nIntCoeffs);
}
int idx = 0;
for (int i = 2; i < P-2; ++i)
{
for (int j = 1; j < Q-i-1; ++j)
{
for (int k = 1; k < R-i-j; ++k)
{
outarray[idx++] = GetMode(i,j,k);
}
}
}
}
int Nektar::StdRegions::StdTetExp::v_GetNedges ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 908 of file StdTetExp.cpp.

{
return 6;
}
int Nektar::StdRegions::StdTetExp::v_GetNfaces ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 913 of file StdTetExp.cpp.

{
return 4;
}
int Nektar::StdRegions::StdTetExp::v_GetNverts ( ) const
protectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 903 of file StdTetExp.cpp.

{
return 4;
}
int Nektar::StdRegions::StdTetExp::v_GetTotalEdgeIntNcoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 986 of file StdTetExp.cpp.

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

{
int P = m_base[0]->GetNumModes()-2;
int Q = m_base[1]->GetNumModes()-2;
int R = m_base[2]->GetNumModes()-2;
return P+Q+4*R;
}
int Nektar::StdRegions::StdTetExp::v_GetTotalFaceIntNcoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1042 of file StdTetExp.cpp.

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

{
int Pi = m_base[0]->GetNumModes() - 2;
int Qi = m_base[1]->GetNumModes() - 2;
int Ri = m_base[2]->GetNumModes() - 2;
return Pi * (2*Qi - Pi - 1) / 2 +
Pi * (2*Ri - Pi - 1) / 2 +
Qi * (2*Ri - Qi - 1);
}
int Nektar::StdRegions::StdTetExp::v_GetVertexMap ( int  localVertexId,
bool  useCoeffPacking = false 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTetExp.

Definition at line 1335 of file StdTetExp.cpp.

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

{
"Mapping not defined for this type of basis");
int localDOF = 0;
if(useCoeffPacking == true) // follow packing of coefficients i.e q,r,p
{
switch(localVertexId)
{
case 0:
{
localDOF = GetMode(0,0,0);
break;
}
case 1:
{
localDOF = GetMode(0,0,1);
break;
}
case 2:
{
localDOF = GetMode(0,1,0);
break;
}
case 3:
{
localDOF = GetMode(1,0,0);
break;
}
default:
{
ASSERTL0(false,"Vertex ID must be between 0 and 3");
break;
}
}
}
else
{
switch(localVertexId)
{
case 0:
{
localDOF = GetMode(0,0,0);
break;
}
case 1:
{
localDOF = GetMode(1,0,0);
break;
}
case 2:
{
localDOF = GetMode(0,1,0);
break;
}
case 3:
{
localDOF = GetMode(0,0,1);
break;
}
default:
{
ASSERTL0(false,"Vertex ID must be between 0 and 3");
break;
}
}
}
return localDOF;
}
void Nektar::StdRegions::StdTetExp::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} (\eta_{1i}) \psi_{pq}^{b} (\eta_{2j}) \psi_{pqr}^{c} (\eta_{3k}) w_i w_j w_k u(\eta_{1,i} \eta_{2,j} \eta_{3,k}) J_{i,j,k}\\ & = & \sum_{i=0}^{nq_0} \psi_p^a(\eta_{1,i}) \sum_{j=0}^{nq_1} \psi_{pq}^b(\eta_{2,j}) \sum_{k=0}^{nq_2} \psi_{pqr}^c u(\eta_{1i},\eta_{2j},\eta_{3k}) J_{i,j,k} \end{array} $

where

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

which can be implemented as
$f_{pqr} (\xi_{3k}) = \sum_{k=0}^{nq_3} \psi_{pqr}^c u(\eta_{1i},\eta_{2j},\eta_{3k}) J_{i,j,k} = {\bf B_3 U} $

$ g_{pq} (\xi_{3k}) = \sum_{j=0}^{nq_1} \psi_{pq}^b (\xi_{2j}) f_{pqr} (\xi_{3k}) = {\bf B_2 F} $

$ (\phi_{pqr}, u)_{\delta} = \sum_{k=0}^{nq_0} \psi_{p}^a (\xi_{3k}) g_{pq} (\xi_{3k}) = {\bf B_1 G} $

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

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTetExp, and Nektar::LocalRegions::TetExp.

Definition at line 506 of file StdTetExp.cpp.

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

Referenced by v_FwdTrans().

{
"Basis[1] is not a general tensor type");
"Basis[2] is not a general tensor type");
if(m_base[0]->Collocation() && m_base[1]->Collocation())
{
MultiplyByQuadratureMetric(inarray,outarray);
}
else
{
}
}
void Nektar::StdRegions::StdTetExp::v_IProductWRTBase_MatOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Definition at line 529 of file StdTetExp.cpp.

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

{
int nq = GetTotPoints();
StdMatrixKey iprodmatkey(eIProductWRTBase,DetShapeType(),*this);
DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
Blas::Dgemv('N',m_ncoeffs,nq,1.0,iprodmat->GetPtr().get(),
m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
}
void Nektar::StdRegions::StdTetExp::v_IProductWRTBase_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual
Parameters
inarrayFunction evaluated at physical collocation points.
outarrayInner product with respect to each basis function over the element.

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTetExp, and Nektar::LocalRegions::TetExp.

Definition at line 548 of file StdTetExp.cpp.

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

Referenced by v_IProductWRTBase().

{
int nquad0 = m_base[0]->GetNumPoints();
int nquad1 = m_base[1]->GetNumPoints();
int nquad2 = m_base[2]->GetNumPoints();
int order0 = m_base[0]->GetNumModes();
int order1 = m_base[1]->GetNumModes();
Array<OneD, NekDouble> tmp (nquad0*nquad1*nquad2);
Array<OneD, NekDouble> wsp (nquad1*nquad2*order0 +
nquad2*order0*(order1+1)/2);
m_base[0]->GetBdata(),
m_base[1]->GetBdata(),
m_base[2]->GetBdata(),
tmp, outarray, wsp, true, true, true);
}
void Nektar::StdRegions::StdTetExp::v_IProductWRTBase_SumFacKernel ( const Array< OneD, const NekDouble > &  base0,
const Array< OneD, const NekDouble > &  base1,
const Array< OneD, const NekDouble > &  base2,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
Array< OneD, NekDouble > &  wsp,
bool  doCheckCollDir0,
bool  doCheckCollDir1,
bool  doCheckCollDir2 
)
protectedvirtual

Implements Nektar::StdRegions::StdExpansion3D.

Definition at line 572 of file StdTetExp.cpp.

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

{
int nquad0 = m_base[0]->GetNumPoints();
int nquad1 = m_base[1]->GetNumPoints();
int nquad2 = m_base[2]->GetNumPoints();
int order0 = m_base[0]->GetNumModes();
int order1 = m_base[1]->GetNumModes();
int order2 = m_base[2]->GetNumModes();
Array<OneD, NekDouble > tmp1 = wsp;
Array<OneD, NekDouble > tmp2 = wsp + nquad1*nquad2*order0;
int i,j, mode,mode1, cnt;
// Inner product with respect to the '0' direction
Blas::Dgemm('T', 'N', nquad1*nquad2, order0, nquad0,
1.0, inarray.get(), nquad0,
base0.get(), nquad0,
0.0, tmp1.get(), nquad1*nquad2);
// Inner product with respect to the '1' direction
for(mode=i=0; i < order0; ++i)
{
Blas::Dgemm('T', 'N', nquad2, order1-i, nquad1,
1.0, tmp1.get()+i*nquad1*nquad2, nquad1,
base1.get()+mode*nquad1, nquad1,
0.0, tmp2.get()+mode*nquad2, nquad2);
mode += order1-i;
}
// fix for modified basis for base singular vertex
{
//base singular vertex and singular edge (1+b)/2
//(1+a)/2 components (makes tmp[nquad2] entry into (1+b)/2)
Blas::Dgemv('T', nquad1, nquad2,
1.0, tmp1.get()+nquad1*nquad2, nquad1,
base1.get()+nquad1, 1,
1.0, tmp2.get()+nquad2, 1);
}
// Inner product with respect to the '2' direction
mode = mode1 = cnt = 0;
for(i = 0; i < order0; ++i)
{
for(j = 0; j < order1-i; ++j, ++cnt)
{
Blas::Dgemv('T', nquad2, order2-i-j,
1.0, base2.get()+mode*nquad2, nquad2,
tmp2.get()+cnt*nquad2, 1,
0.0, outarray.get()+mode1, 1);
mode += order2-i-j;
mode1 += order2-i-j;
}
//increment mode in case order1!=order2
for(j = order1-i; j < order2-i; ++j)
{
mode += order2-i-j;
}
}
// fix for modified basis for top singular vertex component
// Already have evaluated (1+c)/2 (1-b)/2 (1-a)/2
{
// add in (1+c)/2 (1+b)/2 component
outarray[1] += Blas::Ddot(nquad2,base2.get()+nquad2,1,
&tmp2[nquad2],1);
// add in (1+c)/2 (1-b)/2 (1+a)/2 component
outarray[1] += Blas::Ddot(nquad2,base2.get()+nquad2,1,
&tmp2[nquad2*order1],1);
}
}
void Nektar::StdRegions::StdTetExp::v_IProductWRTDerivBase ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTetExp, and Nektar::LocalRegions::TetExp.

Definition at line 658 of file StdTetExp.cpp.

References v_IProductWRTDerivBase_SumFac().

{
}
void Nektar::StdRegions::StdTetExp::v_IProductWRTDerivBase_MatOp ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Definition at line 667 of file StdTetExp.cpp.

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

{
ASSERTL0((dir==0)||(dir==1)||(dir==2),"input dir is out of range");
int nq = GetTotPoints();
MatrixType mtype;
switch (dir)
{
case 0:
break;
case 1:
break;
case 2:
break;
}
StdMatrixKey iprodmatkey(mtype,DetShapeType(),*this);
DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
Blas::Dgemv('N',m_ncoeffs,nq,1.0,iprodmat->GetPtr().get(),
m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
}
void Nektar::StdRegions::StdTetExp::v_IProductWRTDerivBase_SumFac ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual
Parameters
inarrayFunction evaluated at physical collocation points.
outarrayInner product with respect to each basis function over the element.

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTetExp.

Definition at line 704 of file StdTetExp.cpp.

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

Referenced by v_IProductWRTDerivBase().

{
int i;
int nquad0 = m_base[0]->GetNumPoints();
int nquad1 = m_base[1]->GetNumPoints();
int nquad2 = m_base[2]->GetNumPoints();
int nqtot = nquad0*nquad1*nquad2;
int nmodes0 = m_base[0]->GetNumModes();
int nmodes1 = m_base[1]->GetNumModes();
int wspsize = nquad0 + nquad1 + nquad2 + max(nqtot,m_ncoeffs)
+ nquad1*nquad2*nmodes0 + nquad2*nmodes0*(nmodes1+1)/2;
Array<OneD, NekDouble> gfac0(wspsize);
Array<OneD, NekDouble> gfac1(gfac0 + nquad0);
Array<OneD, NekDouble> gfac2(gfac1 + nquad1);
Array<OneD, NekDouble> tmp0 (gfac2 + nquad2);
Array<OneD, NekDouble> wsp(tmp0 + max(nqtot,m_ncoeffs));
const Array<OneD, const NekDouble>& z0 = m_base[0]->GetZ();
const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
const Array<OneD, const NekDouble>& z2 = m_base[2]->GetZ();
// set up geometric factor: (1+z0)/2
for(i = 0; i < nquad0; ++i)
{
gfac0[i] = 0.5*(1+z0[i]);
}
// set up geometric factor: 2/(1-z1)
for(i = 0; i < nquad1; ++i)
{
gfac1[i] = 2.0/(1-z1[i]);
}
// Set up geometric factor: 2/(1-z2)
for(i = 0; i < nquad2; ++i)
{
gfac2[i] = 2.0/(1-z2[i]);
}
// Derivative in first direction is always scaled as follows
for(i = 0; i < nquad1*nquad2; ++i)
{
Vmath::Smul(nquad0,gfac1[i%nquad1],&inarray[0]+i*nquad0,1,&tmp0[0]+i*nquad0,1);
}
for(i = 0; i < nquad2; ++i)
{
Vmath::Smul(nquad0*nquad1,gfac2[i],&tmp0[0]+i*nquad0*nquad1,1,&tmp0[0]+i*nquad0*nquad1,1);
}
switch(dir)
{
case 0:
{
m_base[1]->GetBdata(),
m_base[2]->GetBdata(),
tmp0,outarray,wsp,
false, true, true);
}
break;
case 1:
{
Array<OneD, NekDouble> tmp3(m_ncoeffs);
for(i = 0; i < nquad1*nquad2; ++i)
{
Vmath::Vmul(nquad0,&gfac0[0],1,&tmp0[0]+i*nquad0,1,&tmp0[0]+i*nquad0,1);
}
m_base[1]->GetBdata(),
m_base[2]->GetBdata(),
tmp0,tmp3,wsp,
false, true, true);
for(i = 0; i < nquad2; ++i)
{
Vmath::Smul(nquad0*nquad1,gfac2[i],&inarray[0]+i*nquad0*nquad1,1,&tmp0[0]+i*nquad0*nquad1,1);
}
m_base[1]->GetDbdata(),
m_base[2]->GetBdata(),
tmp0,outarray,wsp,
true, false, true);
Vmath::Vadd(m_ncoeffs,&tmp3[0],1,&outarray[0],1,&outarray[0],1);
}
break;
case 2:
{
Array<OneD, NekDouble> tmp3(m_ncoeffs);
Array<OneD, NekDouble> tmp4(m_ncoeffs);
for(i = 0; i < nquad1; ++i)
{
gfac1[i] = (1+z1[i])/2;
}
for(i = 0; i < nquad1*nquad2; ++i)
{
Vmath::Vmul(nquad0,&gfac0[0],1,&tmp0[0]+i*nquad0,1,&tmp0[0]+i*nquad0,1);
}
m_base[1]->GetBdata(),
m_base[2]->GetBdata(),
tmp0,tmp3,wsp,
false, true, true);
for(i = 0; i < nquad2; ++i)
{
Vmath::Smul(nquad0*nquad1,gfac2[i],&inarray[0]+i*nquad0*nquad1,1,&tmp0[0]+i*nquad0*nquad1,1);
}
for(i = 0; i < nquad1*nquad2; ++i)
{
Vmath::Smul(nquad0,gfac1[i%nquad1],&tmp0[0]+i*nquad0,1,&tmp0[0]+i*nquad0,1);
}
m_base[1]->GetDbdata(),
m_base[2]->GetBdata(),
tmp0,tmp4,wsp,
true, false, true);
m_base[1]->GetBdata(),
m_base[2]->GetDbdata(),
tmp0,outarray,wsp,
true, true, false);
Vmath::Vadd(m_ncoeffs,&tmp3[0],1,&outarray[0],1,&outarray[0],1);
Vmath::Vadd(m_ncoeffs,&tmp4[0],1,&outarray[0],1,&outarray[0],1);
}
break;
default:
{
ASSERTL1(false, "input dir is out of range");
}
break;
}
}
bool Nektar::StdRegions::StdTetExp::v_IsBoundaryInteriorExpansion ( )
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1180 of file StdTetExp.cpp.

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

Referenced by v_GetFaceToElementMap().

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 857 of file StdTetExp.cpp.

References Nektar::NekConstants::kNekZeroTol.

{
if( fabs(xi[2]-1.0) < NekConstants::kNekZeroTol)
{
// Very top point of the tetrahedron
eta[0] = -1.0;
eta[1] = -1.0;
eta[2] = xi[2];
}
else
{
if( fabs(xi[1]-1.0) < NekConstants::kNekZeroTol )
{
// Distant diagonal edge shared by all eta_x
// coordinate planes: the xi_y == -xi_z line
eta[0] = -1.0;
}
else if (fabs(xi[1] + xi[2]) < NekConstants::kNekZeroTol)
{
eta[0] = -1.0;
}
else
{
eta[0] = 2.0*(1.0+xi[0])/(-xi[1]-xi[2]) - 1.0;
}
eta[1] = 2.0*(1.0+xi[1])/(1.0-xi[2]) - 1.0;
eta[2] = xi[2];
}
}
void Nektar::StdRegions::StdTetExp::v_MultiplyByStdQuadratureMetric ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1784 of file StdTetExp.cpp.

References Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauMAlpha2Beta0, Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::StdRegions::StdExpansion::m_base, and Vmath::Vmul().

{
int i, j;
int nquad0 = m_base[0]->GetNumPoints();
int nquad1 = m_base[1]->GetNumPoints();
int nquad2 = m_base[2]->GetNumPoints();
const Array<OneD, const NekDouble>& w0 = m_base[0]->GetW();
const Array<OneD, const NekDouble>& w1 = m_base[1]->GetW();
const Array<OneD, const NekDouble>& w2 = m_base[2]->GetW();
const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
const Array<OneD, const NekDouble>& z2 = m_base[2]->GetZ();
// multiply by integration constants
for(i = 0; i < nquad1*nquad2; ++i)
{
Vmath::Vmul(nquad0,(NekDouble*)&inarray[0]+i*nquad0,1,
w0.get(),1, &outarray[0]+i*nquad0,1);
}
switch(m_base[1]->GetPointsType())
{
// (1,0) Jacobi Inner product.
for(j = 0; j < nquad2; ++j)
{
for(i = 0; i < nquad1; ++i)
{
Blas::Dscal(nquad0,0.5*w1[i], &outarray[0]+i*nquad0+
j*nquad0*nquad1,1);
}
}
break;
default:
for(j = 0; j < nquad2; ++j)
{
for(i = 0; i < nquad1; ++i)
{
Blas::Dscal(nquad0,
0.5*(1-z1[i])*w1[i],
&outarray[0]+i*nquad0 + j*nquad0*nquad1,
1 );
}
}
break;
}
switch(m_base[2]->GetPointsType())
{
// (2,0) Jacobi inner product.
for(i = 0; i < nquad2; ++i)
{
Blas::Dscal(nquad0*nquad1, 0.25*w2[i],
&outarray[0]+i*nquad0*nquad1, 1);
}
break;
default:
for(i = 0; i < nquad2; ++i)
{
Blas::Dscal(nquad0*nquad1,0.25*(1-z2[i])*(1-z2[i])*w2[i],
&outarray[0]+i*nquad0*nquad1,1);
}
break;
}
}
int Nektar::StdRegions::StdTetExp::v_NumBndryCoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 923 of file StdTetExp.cpp.

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

{
"BasisType is not a boundary interior form");
"BasisType is not a boundary interior form");
"BasisType is not a boundary interior form");
int P = m_base[0]->GetNumModes();
int Q = m_base[1]->GetNumModes();
int R = m_base[2]->GetNumModes();
}
int Nektar::StdRegions::StdTetExp::v_NumDGBndryCoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 943 of file StdTetExp.cpp.

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

{
"BasisType is not a boundary interior form");
"BasisType is not a boundary interior form");
"BasisType is not a boundary interior form");
int P = m_base[0]->GetNumModes()-1;
int Q = m_base[1]->GetNumModes()-1;
int R = m_base[2]->GetNumModes()-1;
return (Q+1) + P*(1 + 2*Q - P)/2 // base face
+ (R+1) + P*(1 + 2*R - P)/2 // front face
+ 2*(R+1) + Q*(1 + 2*R - Q); // back two faces
}
void Nektar::StdRegions::StdTetExp::v_PhysDeriv ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_dxi0,
Array< OneD, NekDouble > &  out_dxi1,
Array< OneD, NekDouble > &  out_dxi2 
)
protectedvirtual

Calculate the derivative of the physical points.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TetExp.

Definition at line 106 of file StdTetExp.cpp.

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

Referenced by v_PhysDeriv(), and v_StdPhysDeriv().

{
int Q0 = m_base[0]->GetNumPoints();
int Q1 = m_base[1]->GetNumPoints();
int Q2 = m_base[2]->GetNumPoints();
int Qtot = Q0*Q1*Q2;
// Compute the physical derivative
Array<OneD, NekDouble> out_dEta0(3*Qtot,0.0);
Array<OneD, NekDouble> out_dEta1 = out_dEta0 + Qtot;
Array<OneD, NekDouble> out_dEta2 = out_dEta1 + Qtot;
bool Do_2 = (out_dxi2.num_elements() > 0)? true:false;
bool Do_1 = (out_dxi1.num_elements() > 0)? true:false;
if(Do_2) // Need all local derivatives
{
PhysTensorDeriv(inarray, out_dEta0, out_dEta1, out_dEta2);
}
else if (Do_1) // Need 0 and 1 derivatives
{
PhysTensorDeriv(inarray, out_dEta0, out_dEta1, NullNekDouble1DArray);
}
else // Only need Eta0 derivaitve
{
}
Array<OneD, const NekDouble> eta_0, eta_1, eta_2;
eta_0 = m_base[0]->GetZ();
eta_1 = m_base[1]->GetZ();
eta_2 = m_base[2]->GetZ();
// calculate 2.0/((1-eta_1)(1-eta_2)) Out_dEta0
NekDouble *dEta0 = &out_dEta0[0];
NekDouble fac;
for(int k=0; k< Q2; ++k)
{
for(int j=0; j<Q1; ++j,dEta0+=Q0)
{
Vmath::Smul(Q0,2.0/(1.0-eta_1[j]),dEta0,1,dEta0,1);
}
fac = 1.0/(1.0-eta_2[k]);
Vmath::Smul(Q0*Q1,fac,&out_dEta0[0]+k*Q0*Q1,1,&out_dEta0[0]+k*Q0*Q1,1);
}
if (out_dxi0.num_elements() > 0)
{
// out_dxi1 = 4.0/((1-eta_1)(1-eta_2)) Out_dEta0
Vmath::Smul(Qtot,2.0,out_dEta0,1,out_dxi0,1);
}
if (Do_1||Do_2)
{
Array<OneD, NekDouble> Fac0(Q0);
Vmath::Sadd(Q0,1.0,eta_0,1,Fac0,1);
// calculate 2.0*(1+eta_0)/((1-eta_1)(1-eta_2)) Out_dEta0
for(int k = 0; k < Q1*Q2; ++k)
{
Vmath::Vmul(Q0,&Fac0[0],1,&out_dEta0[0]+k*Q0,1,&out_dEta0[0]+k*Q0,1);
}
// calculate 2/(1.0-eta_2) out_dEta1
for(int k = 0; k < Q2; ++k)
{
Vmath::Smul(Q0*Q1,2.0/(1.0-eta_2[k]),&out_dEta1[0]+k*Q0*Q1,1,
&out_dEta1[0]+k*Q0*Q1,1);
}
if(Do_1)
{
// calculate out_dxi1 = 2.0(1+eta_0)/((1-eta_1)(1-eta_2)) Out_dEta0
// + 2/(1.0-eta_2) out_dEta1
Vmath::Vadd(Qtot,out_dEta0,1,out_dEta1,1,out_dxi1,1);
}
if(Do_2)
{
// calculate (1 + eta_1)/(1 -eta_2)*out_dEta1
NekDouble *dEta1 = &out_dEta1[0];
for(int k=0; k< Q2; ++k)
{
for(int j=0; j<Q1; ++j,dEta1+=Q0)
{
Vmath::Smul(Q0,(1.0+eta_1[j])/2.0,dEta1,1,dEta1,1);
}
}
// calculate out_dxi1 =
// 2.0(1+eta_0)/((1-eta_1)(1-eta_2)) Out_dEta0 +
// (1 + eta_1)/(1 -eta_2)*out_dEta1 + out_dEta2
Vmath::Vadd(Qtot,out_dEta0,1,out_dEta1,1,out_dxi2,1);
Vmath::Vadd(Qtot,out_dEta2,1,out_dxi2 ,1,out_dxi2,1);
}
}
}
void Nektar::StdRegions::StdTetExp::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.

Definition at line 218 of file StdTetExp.cpp.

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

{
switch(dir)
{
case 0:
{
v_PhysDeriv(inarray, outarray, NullNekDouble1DArray,
break;
}
case 1:
{
v_PhysDeriv(inarray, NullNekDouble1DArray, outarray,
break;
}
case 2:
{
break;
}
default:
{
ASSERTL1(false, "input dir is out of range");
}
break;
}
}
void Nektar::StdRegions::StdTetExp::v_StdPhysDeriv ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d0,
Array< OneD, NekDouble > &  out_d1,
Array< OneD, NekDouble > &  out_d2 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 251 of file StdTetExp.cpp.

References v_PhysDeriv().

{
StdTetExp::v_PhysDeriv(inarray, out_d0, out_d1, out_d2);
}
void Nektar::StdRegions::StdTetExp::v_StdPhysDeriv ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 260 of file StdTetExp.cpp.

References v_PhysDeriv().

{
StdTetExp::v_PhysDeriv(dir, inarray, outarray);
}
void Nektar::StdRegions::StdTetExp::v_SVVLaplacianFilter ( Array< OneD, NekDouble > &  array,
const StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1857 of file StdTetExp.cpp.

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

{
//To do : 1) add a test to ensure 0 \leq SvvCutoff \leq 1.
// 2) check if the transfer function needs an analytical
// Fourier transform.
// 3) if it doesn't : find a transfer function that renders
// the if( cutoff_a ...) useless to reduce computational
// cost.
// 4) add SVVDiffCoef to both models!!
int qa = m_base[0]->GetNumPoints();
int qb = m_base[1]->GetNumPoints();
int qc = m_base[2]->GetNumPoints();
int nmodes_a = m_base[0]->GetNumModes();
int nmodes_b = m_base[1]->GetNumModes();
int nmodes_c = m_base[2]->GetNumModes();
// Declare orthogonal basis.
LibUtilities::PointsKey pa(qa,m_base[0]->GetPointsType());
LibUtilities::PointsKey pb(qb,m_base[1]->GetPointsType());
LibUtilities::PointsKey pc(qc,m_base[2]->GetPointsType());
LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A,nmodes_a,pa);
LibUtilities::BasisKey Bb(LibUtilities::eOrtho_B,nmodes_b,pb);
LibUtilities::BasisKey Bc(LibUtilities::eOrtho_C,nmodes_c,pc);
StdTetExp OrthoExp(Ba,Bb,Bc);
Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
int i,j,k,cnt = 0;
//SVV filter paramaters (how much added diffusion relative to physical one
// and fraction of modes from which you start applying this added diffusion)
//
NekDouble SvvDiffCoeff = mkey.GetConstFactor(StdRegions::eFactorSVVDiffCoeff);
NekDouble SVVCutOff = mkey.GetConstFactor(StdRegions::eFactorSVVCutoffRatio);
//Defining the cut of mode
int cutoff_a = (int) (SVVCutOff*nmodes_a);
int cutoff_b = (int) (SVVCutOff*nmodes_b);
int cutoff_c = (int) (SVVCutOff*nmodes_c);
int nmodes = min(min(nmodes_a,nmodes_b),nmodes_c);
NekDouble cutoff = min(min(cutoff_a,cutoff_b),cutoff_c);
NekDouble epsilon = 1;
// project onto physical space.
OrthoExp.FwdTrans(array,orthocoeffs);
//------"New" Version August 22nd '13--------------------
for(i = 0; i < nmodes_a; ++i)
{
for(j = 0; j < nmodes_b-i; ++j)
{
for(k = 0; k < nmodes_c-i-j; ++k)
{
if(i + j + k >= cutoff)
{
orthocoeffs[cnt] *= ((1.0+SvvDiffCoeff)*exp(-(i+j+k-nmodes)*(i+j+k-nmodes)/((NekDouble)((i+j+k-cutoff+epsilon)*(i+j+k-cutoff+epsilon)))));
}
cnt++;
}
}
}
// backward transform to physical space
OrthoExp.BwdTrans(orthocoeffs,array);
}