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

#include <StdTriExp.h>

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

Public Member Functions

 StdTriExp ()
 
 StdTriExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb)
 
 StdTriExp (const StdTriExp &T)
 
 ~StdTriExp ()
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion2D
 StdExpansion2D ()
 
 StdExpansion2D (int numcoeffs, const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb)
 
 StdExpansion2D (const StdExpansion2D &T)
 
virtual ~StdExpansion2D ()
 
void PhysTensorDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d0, Array< OneD, NekDouble > &outarray_d1)
 Calculate the 2D derivative in the local tensor/collapsed coordinate at the physical points. More...
 
NekDouble Integral (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &w0, const Array< OneD, const NekDouble > &w1)
 
void BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
 
void IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
 
- 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
 
boost::shared_ptr< StdExpansionGetLinStdExp (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)
 
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 GetFaceNumModes (const int fid, const Orientation faceOrient, int &numModes0, int &numModes1)
 
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_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 NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray)
 Integrates the specified function over the domain. More...
 
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=NullNekDouble1DArray)
 Calculate the derivative of the physical points. More...
 
virtual void v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Calculate the derivative of the physical points in a given direction. More...
 
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=NullNekDouble1DArray)
 
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)
 Backward tranform for triangular elements. More...
 
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 > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1)
 
virtual void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Transform a given function from physical quadrature space to coefficient space. More...
 
virtual void v_FwdTrans_BndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Calculate the inner product of inarray with respect to the basis B=base0[p]*base1[pq] and put into outarray. More...
 
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 multiplybyweights=true)
 
virtual void v_IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1)
 
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 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_GetEdgeNumPoints (const int i) const
 
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
virtual LibUtilities::BasisType v_GetEdgeBasisType (const int i) const
 
virtual void v_GetCoords (Array< OneD, NekDouble > &coords_x, Array< OneD, NekDouble > &coords_y, Array< OneD, NekDouble > &coords_z)
 
virtual bool v_IsBoundaryInteriorExpansion ()
 
virtual int v_DetCartesianDirOfEdge (const int edge)
 
virtual const
LibUtilities::BasisKey 
v_DetEdgeBasisKey (const int edge) const
 
virtual void v_GetEdgeToElementMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int P=-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_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_SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
virtual void v_ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
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_GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true)
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion2D
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 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::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

Definition at line 50 of file StdTriExp.h.

Constructor & Destructor Documentation

Nektar::StdRegions::StdTriExp::StdTriExp ( )

Definition at line 47 of file StdTriExp.cpp.

48  {
49  }
Nektar::StdRegions::StdTriExp::StdTriExp ( const LibUtilities::BasisKey Ba,
const LibUtilities::BasisKey Bb 
)

Definition at line 52 of file StdTriExp.cpp.

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

54  :
56  Ba.GetNumModes(),
57  Bb.GetNumModes()),
58  2,Ba,Bb),
60  Ba.GetNumModes(),
61  Bb.GetNumModes()),
62  Ba,Bb)
63  {
64  ASSERTL0(Ba.GetNumModes() <= Bb.GetNumModes(),
65  "order in 'a' direction is higher than order "
66  "in 'b' direction");
67  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
int getNumberOfCoefficients(int Na, int Nb)
Definition: ShapeType.hpp:111
StdExpansion()
Default Constructor.
Nektar::StdRegions::StdTriExp::StdTriExp ( const StdTriExp T)

Definition at line 69 of file StdTriExp.cpp.

69  :
70  StdExpansion(T),
72  {
73  }
StdExpansion()
Default Constructor.
Nektar::StdRegions::StdTriExp::~StdTriExp ( )

Definition at line 75 of file StdTriExp.cpp.

76  {
77  }

Member Function Documentation

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

Backward tranform for triangular elements.

Note
'q' (base[1]) runs fastest in this element.

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp.

Definition at line 246 of file StdTriExp.cpp.

References v_BwdTrans_SumFac().

249  {
250  v_BwdTrans_SumFac(inarray,outarray);
251  }
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdTriExp.cpp:254
void Nektar::StdRegions::StdTriExp::v_BwdTrans_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::NodalTriExp, and Nektar::StdRegions::StdNodalTriExp.

Definition at line 254 of file StdTriExp.cpp.

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

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

257  {
258  Array<OneD, NekDouble> wsp(m_base[0]->GetNumPoints()*
259  m_base[1]->GetNumModes());
260 
261  BwdTrans_SumFacKernel(m_base[0]->GetBdata(),
262  m_base[1]->GetBdata(),
263  inarray,outarray,wsp);
264  }
void BwdTrans_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
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::StdTriExp::v_BwdTrans_SumFacKernel ( const Array< OneD, const NekDouble > &  base0,
const Array< OneD, const NekDouble > &  base1,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
Array< OneD, NekDouble > &  wsp,
bool  doCheckCollDir0,
bool  doCheckCollDir1 
)
protectedvirtual

Implements Nektar::StdRegions::StdExpansion2D.

Definition at line 266 of file StdTriExp.cpp.

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

274  {
275  int i;
276  int mode;
277  int nquad0 = m_base[0]->GetNumPoints();
278  int nquad1 = m_base[1]->GetNumPoints();
279  int nmodes0 = m_base[0]->GetNumModes();
280  int nmodes1 = m_base[1]->GetNumModes();
281 
282  ASSERTL1(wsp.num_elements() >= nquad0*nmodes1,
283  "Workspace size is not sufficient");
286  "Basis[1] is not of general tensor type");
287 
288  for (i = mode = 0; i < nmodes0; ++i)
289  {
290  Blas::Dgemv('N', nquad1,nmodes1-i,1.0,base1.get()+mode*nquad1,
291  nquad1,&inarray[0]+mode,1,0.0,&wsp[0]+i*nquad1,1);
292  mode += nmodes1-i;
293  }
294 
295  // fix for modified basis by splitting top vertex mode
297  {
298  Blas::Daxpy(nquad1,inarray[1],base1.get()+nquad1,1,
299  &wsp[0]+nquad1,1);
300  }
301 
302  Blas::Dgemm('N','T', nquad0,nquad1,nmodes0,1.0, base0.get(),nquad0,
303  &wsp[0], nquad1,0.0, &outarray[0], nquad0);
304  }
Principle Modified Functions .
Definition: BasisType.h:49
Principle Orthogonal Functions .
Definition: BasisType.h:47
Principle Modified Functions .
Definition: BasisType.h:50
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:250
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
Array< OneD, LibUtilities::BasisSharedPtr > m_base
int Nektar::StdRegions::StdTriExp::v_CalcNumberOfCoefficients ( const std::vector< unsigned int > &  nummodes,
int &  modes_offset 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 798 of file StdTriExp.cpp.

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

801  {
803  nummodes[modes_offset],
804  nummodes[modes_offset+1]);
805  modes_offset += 2;
806 
807  return nmodes;
808  }
int getNumberOfCoefficients(int Na, int Nb)
Definition: ShapeType.hpp:111
DNekMatSharedPtr Nektar::StdRegions::StdTriExp::v_CreateStdMatrix ( const StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::NodalTriExp, Nektar::LocalRegions::TriExp, and Nektar::StdRegions::StdNodalTriExp.

Definition at line 1329 of file StdTriExp.cpp.

References v_GenMatrix().

1330  {
1331  return v_GenMatrix(mkey);
1332  }
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey)
Definition: StdTriExp.cpp:1252
int Nektar::StdRegions::StdTriExp::v_DetCartesianDirOfEdge ( const int  edge)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 851 of file StdTriExp.cpp.

References ASSERTL2.

852  {
853  ASSERTL2(edge >= 0 && edge <= 2, "edge id is out of range");
854 
855  return edge == 0 ? 0 : 1;
856  }
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:250
const LibUtilities::BasisKey Nektar::StdRegions::StdTriExp::v_DetEdgeBasisKey ( const int  edge) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 858 of file StdTriExp.cpp.

References ASSERTL0, ASSERTL2, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::StdRegions::StdExpansion::GetBasis(), Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::GetNumPoints(), Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::StdRegions::StdExpansion::m_base, and Nektar::LibUtilities::NullBasisKey().

860  {
861  ASSERTL2(i >= 0 && i <= 2, "edge id is out of range");
862 
863  if (i == 0)
864  {
865  return GetBasis(0)->GetBasisKey();
866  }
867  else
868  {
869  switch(m_base[1]->GetBasisType())
870  {
872  {
873  switch(m_base[1]->GetPointsType())
874  {
876  {
877  LibUtilities::PointsKey pkey(
878  m_base[1]->GetBasisKey().GetPointsKey().
879  GetNumPoints()+1,
881  return LibUtilities::BasisKey(
883  m_base[1]->GetNumModes(),pkey);
884  break;
885  }
886 
887  default:
888  ASSERTL0(false,"unexpected points distribution");
889  break;
890  }
891  }
892  default:
893  ASSERTL0(false,"Information not available to set edge key");
894  break;
895  }
896  }
898  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
const LibUtilities::BasisSharedPtr & GetBasis(int dir) const
This function gets the shared point to basis in the dir direction.
Definition: StdExpansion.h:118
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:229
Principle Modified Functions .
Definition: BasisType.h:49
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:59
Principle Modified Functions .
Definition: BasisType.h:50
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:250
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.
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
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:52
LibUtilities::ShapeType Nektar::StdRegions::StdTriExp::v_DetShapeType ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 745 of file StdTriExp.cpp.

References Nektar::LibUtilities::eTriangle.

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp.

Definition at line 686 of file StdTriExp.cpp.

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

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

688  {
689  int i,m;
690  int nquad0 = m_base[0]->GetNumPoints();
691  int nquad1 = m_base[1]->GetNumPoints();
692  int order0 = m_base[0]->GetNumModes();
693  int order1 = m_base[1]->GetNumModes();
694  int mode0 = 0;
695  Array<OneD, const NekDouble> base0 = m_base[0]->GetBdata();
696  Array<OneD, const NekDouble> base1 = m_base[1]->GetBdata();
697 
698  ASSERTL2(mode <= m_ncoeffs,
699  "calling argument mode is larger than "
700  "total expansion order");
701 
702  m = order1;
703  for (i = 0; i < order0; ++i, m+=order1-i)
704  {
705  if (m > mode)
706  {
707  mode0 = i;
708  break;
709  }
710  }
711 
712  // deal with top vertex mode in modified basis
713  if (mode == 1 &&
715  {
716  Vmath::Fill(nquad0*nquad1 , 1.0, outarray, 1);
717  }
718  else
719  {
720  for (i = 0; i < nquad1; ++i)
721  {
722  Vmath::Vcopy(nquad0,(NekDouble *)(base0.get()+mode0*nquad0),
723  1,&outarray[0]+i*nquad0,1);
724  }
725  }
726 
727  for(i = 0; i < nquad0; ++i)
728  {
729  Vmath::Vmul(nquad1,(NekDouble *)(base1.get() + mode*nquad1),
730  1,&outarray[0]+i,nquad0,&outarray[0]+i,nquad0);
731  }
732  }
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:46
Principle Modified Functions .
Definition: BasisType.h:49
double NekDouble
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:250
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
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:183
void Nektar::StdRegions::StdTriExp::v_FwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Transform a given function from physical quadrature space to coefficient space.

See also
StdExpansion::FwdTrans

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::NodalTriExp, Nektar::LocalRegions::TriExp, and Nektar::StdRegions::StdNodalTriExp.

Definition at line 306 of file StdTriExp.cpp.

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

309  {
310  v_IProductWRTBase(inarray,outarray);
311 
312  // get Mass matrix inverse
313  StdMatrixKey masskey(eInvMass,DetShapeType(),*this);
314  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
315 
316  // copy inarray in case inarray == outarray
317  NekVector<NekDouble> in(m_ncoeffs,outarray,eCopy);
318  NekVector<NekDouble> out(m_ncoeffs,outarray,eWrapper);
319 
320  out = (*matsys)*in;
321  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:470
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Calculate the inner product of inarray with respect to the basis B=base0[p]*base1[pq] and put into ou...
Definition: StdTriExp.cpp:460
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:706
void Nektar::StdRegions::StdTriExp::v_FwdTrans_BndConstrained ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp.

Definition at line 324 of file StdTriExp.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::StdRegions::eForwards, Nektar::StdRegions::eMass, Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::StdRegions::StdExpansion::m_stdStaticCondMatrixManager, Nektar::StdRegions::StdExpansion::MassMatrixOp(), sign, v_GetEdgeToElementMap(), v_GetInteriorMap(), v_IProductWRTBase(), v_NumBndryCoeffs(), and Vmath::Vsub().

327  {
328  int i,j;
329  int npoints[2] = {m_base[0]->GetNumPoints(),
330  m_base[1]->GetNumPoints()};
331  int nmodes[2] = {m_base[0]->GetNumModes(),
332  m_base[1]->GetNumModes()};
333 
334  fill(outarray.get(), outarray.get()+m_ncoeffs, 0.0 );
335 
336  Array<OneD, NekDouble> physEdge[3];
337  Array<OneD, NekDouble> coeffEdge[3];
338  for(i = 0; i < 3; i++)
339  {
340  physEdge[i] = Array<OneD, NekDouble>(npoints[i!=0]);
341  coeffEdge[i] = Array<OneD, NekDouble>(nmodes[i!=0]);
342  }
343 
344  for(i = 0; i < npoints[0]; i++)
345  {
346  physEdge[0][i] = inarray[i];
347  }
348 
349  for(i = 0; i < npoints[1]; i++)
350  {
351  physEdge[1][i] = inarray[npoints[0]-1+i*npoints[0]];
352  physEdge[2][i] = inarray[(npoints[1]-1)*npoints[0]-i*npoints[0]];
353  }
354 
355  StdSegExpSharedPtr segexp[2] = {
357  m_base[0]->GetBasisKey()),
359  m_base[1]->GetBasisKey())
360  };
361 
362  Array<OneD, unsigned int> mapArray;
363  Array<OneD, int> signArray;
364  NekDouble sign;
365 
366  for (i = 0; i < 3; i++)
367  {
368  //segexp[i!=0]->v_FwdTrans_BndConstrained(physEdge[i],coeffEdge[i]);
369  segexp[i!=0]->FwdTrans_BndConstrained(physEdge[i],coeffEdge[i]);
370 
371  v_GetEdgeToElementMap(i,eForwards,mapArray,signArray);
372  for (j = 0; j < nmodes[i != 0]; j++)
373  {
374  sign = (NekDouble) signArray[j];
375  outarray[ mapArray[j] ] = sign * coeffEdge[i][j];
376  }
377  }
378 
379  Array<OneD, NekDouble> tmp0(m_ncoeffs);
380  Array<OneD, NekDouble> tmp1(m_ncoeffs);
381 
382  StdMatrixKey masskey(eMass,DetShapeType(),*this);
383  MassMatrixOp(outarray,tmp0,masskey);
384  v_IProductWRTBase(inarray,tmp1);
385 
386  Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
387 
388  // get Mass matrix inverse (only of interior DOF)
389  // use block (1,1) of the static condensed system
390  // note: this block alreay contains the inverse matrix
391  DNekMatSharedPtr matsys =
392  (m_stdStaticCondMatrixManager[masskey])->GetBlock(1,1);
393 
394  int nBoundaryDofs = v_NumBndryCoeffs();
395  int nInteriorDofs = m_ncoeffs - nBoundaryDofs;
396 
397  Array<OneD, NekDouble> rhs (nInteriorDofs);
398  Array<OneD, NekDouble> result(nInteriorDofs);
399 
400  v_GetInteriorMap(mapArray);
401 
402  for (i = 0; i < nInteriorDofs; i++)
403  {
404  rhs[i] = tmp1[ mapArray[i] ];
405  }
406 
407  Blas::Dgemv('N',nInteriorDofs,nInteriorDofs,
408  1.0,&(matsys->GetPtr())[0],nInteriorDofs,
409  rhs.get(),1,
410  0.0,result.get(),1);
411 
412  for (i = 0; i < nInteriorDofs; i++)
413  {
414  outarray[ mapArray[i] ] = result[i];
415  }
416  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:470
virtual void v_GetEdgeToElementMap(const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int P=-1)
Definition: StdTriExp.cpp:906
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:976
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:27
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Calculate the inner product of inarray with respect to the basis B=base0[p]*base1[pq] and put into ou...
Definition: StdTriExp.cpp:460
LibUtilities::NekManager< StdMatrixKey, DNekBlkMat, StdMatrixKey::opLess > m_stdStaticCondMatrixManager
boost::shared_ptr< StdSegExp > StdSegExpSharedPtr
Definition: StdSegExp.h:47
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
virtual void v_GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdTriExp.cpp:1185
virtual int v_NumBndryCoeffs() const
Definition: StdTriExp.cpp:750
double NekDouble
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:343
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdTriExp::v_GeneralMatrixOp_MatOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdMatrixKey mkey 
)
protectedvirtual

Reimplemented in Nektar::LocalRegions::TriExp.

Definition at line 1533 of file StdTriExp.cpp.

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

1537  {
1539 
1540  if(inarray.get() == outarray.get())
1541  {
1542  Array<OneD,NekDouble> tmp(m_ncoeffs);
1543  Vmath::Vcopy(m_ncoeffs,inarray.get(),1,tmp.get(),1);
1544 
1545  Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, 1.0, mat->GetPtr().get(),
1546  m_ncoeffs, tmp.get(), 1, 0.0, outarray.get(), 1);
1547  }
1548  else
1549  {
1550  Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, 1.0, mat->GetPtr().get(),
1551  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
1552  }
1553  }
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:1061
DNekMatSharedPtr Nektar::StdRegions::StdTriExp::v_GenMatrix ( const StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp, Nektar::StdRegions::StdNodalTriExp, and Nektar::LocalRegions::NodalTriExp.

Definition at line 1252 of file StdTriExp.cpp.

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

Referenced by v_CreateStdMatrix().

1253  {
1254 
1255  MatrixType mtype = mkey.GetMatrixType();
1256 
1257  DNekMatSharedPtr Mat;
1258 
1259  switch(mtype)
1260  {
1262  {
1263  int nq0, nq1, nq;
1264 
1265  nq0 = m_base[0]->GetNumPoints();
1266  nq1 = m_base[1]->GetNumPoints();
1267 
1268  // take definition from key
1269  if(mkey.ConstFactorExists(eFactorConst))
1270  {
1271  nq = (int) mkey.GetConstFactor(eFactorConst);
1272  }
1273  else
1274  {
1275  nq = max(nq0,nq1);
1276  }
1277 
1278  int neq = LibUtilities::StdTriData::
1280  Array<OneD, Array<OneD, NekDouble> > coords(neq);
1281  Array<OneD, NekDouble> coll (2);
1282  Array<OneD, DNekMatSharedPtr> I (2);
1283  Array<OneD, NekDouble> tmp (nq0);
1284 
1285 
1286  Mat = MemoryManager<DNekMat>::AllocateSharedPtr(neq,nq0*nq1);
1287  int cnt = 0;
1288 
1289  for(int i = 0; i < nq; ++i)
1290  {
1291  for(int j = 0; j < nq-i; ++j,++cnt)
1292  {
1293  coords[cnt] = Array<OneD, NekDouble>(2);
1294  coords[cnt][0] = -1.0 + 2*j/(NekDouble)(nq-1);
1295  coords[cnt][1] = -1.0 + 2*i/(NekDouble)(nq-1);
1296  }
1297  }
1298 
1299  for(int i = 0; i < neq; ++i)
1300  {
1301  LocCoordToLocCollapsed(coords[i],coll);
1302 
1303  I[0] = m_base[0]->GetI(coll);
1304  I[1] = m_base[1]->GetI(coll+1);
1305 
1306  // interpolate first coordinate direction
1307  for (int j = 0; j < nq1; ++j)
1308  {
1309  NekDouble fac = (I[1]->GetPtr())[j];
1310  Vmath::Smul(nq0, fac, I[0]->GetPtr(), 1, tmp, 1);
1311 
1312  Vmath::Vcopy(nq0, &tmp[0], 1,
1313  Mat->GetRawPtr() + j*nq0*neq + i, neq);
1314  }
1315 
1316  }
1317  break;
1318  }
1319  default:
1320  {
1322  break;
1323  }
1324  }
1325 
1326  return Mat;
1327  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
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:213
int getNumberOfCoefficients(int Na, int Nb)
Definition: ShapeType.hpp:111
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
double NekDouble
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
void Nektar::StdRegions::StdTriExp::v_GetBoundaryMap ( Array< OneD, unsigned int > &  outarray)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp.

Definition at line 1215 of file StdTriExp.cpp.

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

1216  {
1219  "Expansion not of expected type");
1220  int i;
1221  int cnt;
1222  int nummodes0, nummodes1;
1223  int value;
1224 
1225  if (outarray.num_elements()!=NumBndryCoeffs())
1226  {
1227  outarray = Array<OneD, unsigned int>(NumBndryCoeffs());
1228  }
1229 
1230  nummodes0 = m_base[0]->GetNumModes();
1231  nummodes1 = m_base[1]->GetNumModes();
1232 
1233  value = 2*nummodes1-1;
1234  for(i = 0; i < value; i++)
1235  {
1236  outarray[i]=i;
1237  }
1238  cnt = value;
1239 
1240  for(i = 0; i < nummodes0-2; i++)
1241  {
1242  outarray[cnt++]=value;
1243  value += nummodes1-2-i;
1244  }
1245  }
Principle Modified Functions .
Definition: BasisType.h:49
Principle Modified Functions .
Definition: BasisType.h:50
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:228
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdTriExp::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::NodalTriExp, and Nektar::LocalRegions::TriExp.

Definition at line 825 of file StdTriExp.cpp.

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

828  {
829  Array<OneD, const NekDouble> z0 = m_base[0]->GetZ();
830  Array<OneD, const NekDouble> z1 = m_base[1]->GetZ();
831  int nq0 = GetNumPoints(0);
832  int nq1 = GetNumPoints(1);
833  int i,j;
834 
835  for(i = 0; i < nq1; ++i)
836  {
837  for(j = 0; j < nq0; ++j)
838  {
839  coords_0[i*nq0+j] = (1+z0[j])*(1-z1[i])/2.0 - 1.0;
840  }
841  Vmath::Fill(nq0,z1[i],&coords_1[0] + i*nq0,1);
842  }
843  }
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:46
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::StdTriExp::v_GetEdgeBasisType ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 810 of file StdTriExp.cpp.

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

811  {
812  ASSERTL2(i >= 0 && i <= 2, "edge id is out of range");
813 
814  if (i == 0)
815  {
816  return GetBasisType(0);
817  }
818  else
819  {
820  return GetBasisType(1);
821  }
822  }
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:250
void Nektar::StdRegions::StdTriExp::v_GetEdgeInteriorMap ( const int  eid,
const Orientation  edgeOrient,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp.

Definition at line 1099 of file StdTriExp.cpp.

References ASSERTL0, Nektar::StdRegions::eBackwards, Nektar::StdRegions::eForwards, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::StdRegions::StdExpansion::GetEdgeBasisType(), Nektar::StdRegions::StdExpansion::GetEdgeNcoeffs(), and Nektar::StdRegions::StdExpansion::m_base.

1104  {
1107  "Mapping not defined for this type of basis");
1108  int i;
1109  const int nummodes1 = m_base[1]->GetNumModes();
1110  const int nEdgeIntCoeffs = GetEdgeNcoeffs(eid)-2;
1111 
1112  if(maparray.num_elements() != nEdgeIntCoeffs)
1113  {
1114  maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
1115  }
1116 
1117  if(signarray.num_elements() != nEdgeIntCoeffs)
1118  {
1119  signarray = Array<OneD, int>(nEdgeIntCoeffs,1);
1120  }
1121  else
1122  {
1123  fill( signarray.get() , signarray.get()+nEdgeIntCoeffs, 1 );
1124  }
1125 
1126  switch(eid)
1127  {
1128  case 0:
1129  {
1130  int cnt = 2*nummodes1 - 1;
1131  for(i = 0; i < nEdgeIntCoeffs; cnt+=nummodes1-2-i, ++i)
1132  {
1133  maparray[i] = cnt;
1134  }
1135 
1136  if(edgeOrient==eBackwards)
1137  {
1138  for(i = 1; i < nEdgeIntCoeffs; i+=2)
1139  {
1140  signarray[i] = -1;
1141  }
1142  }
1143  break;
1144  }
1145  case 1:
1146  {
1147  for(i = 0; i < nEdgeIntCoeffs; i++)
1148  {
1149  maparray[i] = nummodes1+1+i;
1150  }
1151 
1152  if(edgeOrient==eBackwards)
1153  {
1154  for(i = 1; i < nEdgeIntCoeffs; i+=2)
1155  {
1156  signarray[i] = -1;
1157  }
1158  }
1159  break;
1160  }
1161  case 2:
1162  {
1163  for(i = 0; i < nEdgeIntCoeffs; i++)
1164  {
1165  maparray[i] = 2+i;
1166  }
1167 
1168  if(edgeOrient==eForwards)
1169  {
1170  for(i = 1; i < nEdgeIntCoeffs; i+=2)
1171  {
1172  signarray[i] = -1;
1173  }
1174  }
1175  break;
1176  }
1177  default:
1178  {
1179  ASSERTL0(false,"eid must be between 0 and 2");
1180  break;
1181  }
1182  }
1183  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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 GetEdgeBasisType(const int i) const
This function returns the type of expansion basis on the i-th edge.
Definition: StdExpansion.h:413
Principle Modified Functions .
Definition: BasisType.h:50
Array< OneD, LibUtilities::BasisSharedPtr > m_base
int Nektar::StdRegions::StdTriExp::v_GetEdgeNcoeffs ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 770 of file StdTriExp.cpp.

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

771  {
772  ASSERTL2(i >= 0 && i <= 2, "edge id is out of range");
773 
774  if (i == 0)
775  {
776  return GetBasisNumModes(0);
777  }
778  else
779  {
780  return GetBasisNumModes(1);
781  }
782  }
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:250
int Nektar::StdRegions::StdTriExp::v_GetEdgeNumPoints ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 784 of file StdTriExp.cpp.

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

785  {
786  ASSERTL2((i >= 0)&&(i <= 2),"edge id is out of range");
787 
788  if (i == 0)
789  {
790  return GetNumPoints(0);
791  }
792  else
793  {
794  return GetNumPoints(1);
795  }
796  }
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:229
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:250
void Nektar::StdRegions::StdTriExp::v_GetEdgeToElementMap ( const int  eid,
const Orientation  edgeOrient,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray,
int  P = -1 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp.

Definition at line 906 of file StdTriExp.cpp.

References ASSERTL0, ASSERTL1, Nektar::StdRegions::eBackwards, Nektar::StdRegions::eForwards, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::StdRegions::StdExpansion::GetEdgeBasisType(), Nektar::StdRegions::StdExpansion::m_base, and class_topology::P.

Referenced by v_FwdTrans_BndConstrained().

912  {
915  "Mapping not defined for this type of basis");
916 
917  int i;
918  int numModes=0;
919  int order0 = m_base[0]->GetNumModes();
920  int order1 = m_base[1]->GetNumModes();
921 
922  switch (eid)
923  {
924  case 0:
925  numModes = order0;
926  break;
927  case 1:
928  case 2:
929  numModes = order1;
930  break;
931  }
932 
933  bool checkForZeroedModes = false;
934  if (P == -1)
935  {
936  P = numModes;
937  }
938  else if(P != numModes)
939  {
940  checkForZeroedModes = true;
941  }
942 
943 
944  if(maparray.num_elements() != P)
945  {
946  maparray = Array<OneD, unsigned int>(P);
947  }
948 
949  if(signarray.num_elements() != P)
950  {
951  signarray = Array<OneD, int>(P,1);
952  }
953  else
954  {
955  fill(signarray.get() , signarray.get()+P, 1);
956  }
957 
958  switch(eid)
959  {
960  case 0:
961  {
962  int cnt = 0;
963  for(i = 0; i < P; cnt+=order1-i, ++i)
964  {
965  maparray[i] = cnt;
966  }
967 
968  if(edgeOrient==eBackwards)
969  {
970  swap( maparray[0] , maparray[1] );
971 
972  for(i = 3; i < P; i+=2)
973  {
974  signarray[i] = -1;
975  }
976  }
977  break;
978  }
979  case 1:
980  {
981  maparray[0] = order1;
982  maparray[1] = 1;
983  for(i = 2; i < P; i++)
984  {
985  maparray[i] = order1-1+i;
986  }
987 
988  if(edgeOrient==eBackwards)
989  {
990  swap( maparray[0] , maparray[1] );
991 
992  for(i = 3; i < P; i+=2)
993  {
994  signarray[i] = -1;
995  }
996  }
997  break;
998  }
999  case 2:
1000  {
1001  for(i = 0; i < P; i++)
1002  {
1003  maparray[i] = i;
1004  }
1005 
1006  if(edgeOrient==eForwards)
1007  {
1008  swap( maparray[0] , maparray[1] );
1009 
1010  for(i = 3; i < P; i+=2)
1011  {
1012  signarray[i] = -1;
1013  }
1014  }
1015  break;
1016  }
1017  default:
1018  ASSERTL0(false,"eid must be between 0 and 2");
1019  break;
1020  }
1021 
1022 
1023  if (checkForZeroedModes)
1024  {
1025  // Zero signmap and set maparray to zero if
1026  // elemental modes are not as large as face modes
1027  for (int j = numModes; j < P; j++)
1028  {
1029  signarray[j] = 0.0;
1030  maparray[j] = maparray[0];
1031  }
1032  }
1033  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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
Principle Modified Functions .
Definition: BasisType.h:50
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdTriExp::v_GetInteriorMap ( Array< OneD, unsigned int > &  outarray)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp.

Definition at line 1185 of file StdTriExp.cpp.

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

Referenced by v_FwdTrans_BndConstrained().

1186  {
1189  "Expansion not of a proper type");
1190 
1191  int i,j;
1192  int cnt = 0;
1193  int nummodes0, nummodes1;
1194  int startvalue;
1195  if(outarray.num_elements()!=GetNcoeffs()-NumBndryCoeffs())
1196  {
1197  outarray = Array<OneD, unsigned int>(GetNcoeffs()-NumBndryCoeffs());
1198  }
1199 
1200  nummodes0 = m_base[0]->GetNumModes();
1201  nummodes1 = m_base[1]->GetNumModes();
1202 
1203  startvalue = 2*nummodes1;
1204 
1205  for(i = 0; i < nummodes0-2; i++)
1206  {
1207  for(j = 0; j < nummodes1-3-i; j++)
1208  {
1209  outarray[cnt++]=startvalue+j;
1210  }
1211  startvalue+=nummodes1-2-i;
1212  }
1213  }
Principle Modified Functions .
Definition: BasisType.h:49
Principle Modified Functions .
Definition: BasisType.h:50
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:131
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
Array< OneD, LibUtilities::BasisSharedPtr > m_base
int Nektar::StdRegions::StdTriExp::v_GetNedges ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 740 of file StdTriExp.cpp.

741  {
742  return 3;
743  }
int Nektar::StdRegions::StdTriExp::v_GetNverts ( ) const
protectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 735 of file StdTriExp.cpp.

736  {
737  return 3;
738  }
void Nektar::StdRegions::StdTriExp::v_GetSimplexEquiSpacedConnectivity ( Array< OneD, int > &  conn,
bool  standard = true 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1604 of file StdTriExp.cpp.

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

1607  {
1608  int np1 = m_base[0]->GetNumPoints();
1609  int np2 = m_base[1]->GetNumPoints();
1610  int np = max(np1,np2);
1611 
1612  conn = Array<OneD, int>(3*(np-1)*(np-1));
1613 
1614  int row = 0;
1615  int rowp1 = 0;
1616  int cnt = 0;
1617  for(int i = 0; i < np-1; ++i)
1618  {
1619  rowp1 += np-i;
1620  for(int j = 0; j < np-i-2; ++j)
1621  {
1622  conn[cnt++] = row +j;
1623  conn[cnt++] = row +j+1;
1624  conn[cnt++] = rowp1 +j;
1625 
1626  conn[cnt++] = rowp1 +j+1;
1627  conn[cnt++] = rowp1 +j;
1628  conn[cnt++] = row +j+1;
1629  }
1630 
1631  conn[cnt++] = row +np-i-2;
1632  conn[cnt++] = row +np-i-1;
1633  conn[cnt++] = rowp1+np-i-2;
1634 
1635  row += np-i;
1636  }
1637  }
Array< OneD, LibUtilities::BasisSharedPtr > m_base
int Nektar::StdRegions::StdTriExp::v_GetVertexMap ( int  localVertexId,
bool  useCoeffPacking = false 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp.

Definition at line 1035 of file StdTriExp.cpp.

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

1036  {
1037  ASSERTL0(
1038  GetEdgeBasisType(localVertexId) == LibUtilities::eModified_A ||
1039  GetEdgeBasisType(localVertexId) == LibUtilities::eModified_B,
1040  "Mapping not defined for this type of basis");
1041 
1042  int localDOF = 0;
1043  if(useCoeffPacking == true)
1044  {
1045  switch(localVertexId)
1046  {
1047  case 0:
1048  {
1049  localDOF = 0;
1050  break;
1051  }
1052  case 1:
1053  {
1054  localDOF = 1;
1055  break;
1056  }
1057  case 2:
1058  {
1059  localDOF = m_base[1]->GetNumModes();
1060  break;
1061  }
1062  default:
1063  {
1064  ASSERTL0(false,"eid must be between 0 and 2");
1065  break;
1066  }
1067  }
1068  }
1069  else // follow book format for vertex indexing.
1070  {
1071  switch(localVertexId)
1072  {
1073  case 0:
1074  {
1075  localDOF = 0;
1076  break;
1077  }
1078  case 1:
1079  {
1080  localDOF = m_base[1]->GetNumModes();
1081  break;
1082  }
1083  case 2:
1084  {
1085  localDOF = 1;
1086  break;
1087  }
1088  default:
1089  {
1090  ASSERTL0(false,"eid must be between 0 and 2");
1091  break;
1092  }
1093  }
1094  }
1095 
1096  return localDOF;
1097  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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
Principle Modified Functions .
Definition: BasisType.h:50
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdTriExp::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::NodalTriExp, Nektar::LocalRegions::TriExp, and Nektar::StdRegions::StdNodalTriExp.

Definition at line 1375 of file StdTriExp.cpp.

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

1379  {
1380  StdTriExp::v_HelmholtzMatrixOp_MatFree(inarray,outarray,mkey);
1381  }
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
NekDouble Nektar::StdRegions::StdTriExp::v_Integral ( const Array< OneD, const NekDouble > &  inarray)
protectedvirtual

Integrates the specified function over the domain.

See also
StdRegions::StdExpansion::Integral.

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::NodalTriExp, and Nektar::LocalRegions::TriExp.

Definition at line 82 of file StdTriExp.cpp.

References Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::StdRegions::StdExpansion2D::Integral(), Nektar::StdRegions::StdExpansion::m_base, and Vmath::Smul().

84  {
85  int i;
86  int nquad1 = m_base[1]->GetNumPoints();
87  Array<OneD, NekDouble> w1_tmp(nquad1);
88 
89  Array<OneD, const NekDouble> w0 = m_base[0]->GetW();
90  Array<OneD, const NekDouble> w1 = m_base[1]->GetW();
91  Array<OneD, const NekDouble> z1 = m_base[1]->GetZ();
92 
93  switch(m_base[1]->GetPointsType())
94  {
95  case LibUtilities::eGaussRadauMAlpha1Beta0: // (0,1) Jacobi Inner product
96  {
97  Vmath::Smul(nquad1, 0.5, w1, 1, w1_tmp,1);
98  break;
99  }
100  default:
101  {
102  // include jacobian factor on whatever coordinates are defined.
103  for(i = 0; i < nquad1; ++i)
104  {
105  w1_tmp[i] = 0.5*(1-z1[i])*w1[i];
106  }
107  break;
108  }
109  }
110 
111  return StdExpansion2D::Integral(inarray,w0,w1_tmp);
112  }
NekDouble Integral(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &w0, const Array< OneD, const NekDouble > &w1)
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:59
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:213
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::StdTriExp::v_IProductWRTBase ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Calculate the inner product of inarray with respect to the basis B=base0[p]*base1[pq] and put into outarray.

$ \begin{array}{rcl} I_{pq} = (\phi^A_q \phi^B_{pq}, u) &=& \sum_{i=0}^{nq_0}\sum_{j=0}^{nq_1} \phi^A_p(\eta_{0,i})\phi^B_{pq}(\eta_{1,j}) w^0_i w^1_j u(\xi_{0,i} \xi_{1,j}) \\ & = & \sum_{i=0}^{nq_0} \phi^A_p(\eta_{0,i}) \sum_{j=0}^{nq_1} \phi^B_{pq}(\eta_{1,j}) \tilde{u}_{i,j} \end{array} $

where

$ \tilde{u}_{i,j} = w^0_i w^1_j u(\xi_{0,i},\xi_{1,j}) $

which can be implemented as

$ f_{pj} = \sum_{i=0}^{nq_0} \phi^A_p(\eta_{0,i}) \tilde{u}_{i,j} \rightarrow {\bf B_1 U} $ $ I_{pq} = \sum_{j=0}^{nq_1} \phi^B_{pq}(\eta_{1,j}) f_{pj} \rightarrow {\bf B_2[p*skip] f[skip]} $

Recall: $ \eta_{1} = \frac{2(1+\xi_1)}{(1-\xi_2)}-1, \, \eta_2 = \xi_2$

Note: For the orthgonality of this expansion to be realised the 'q' ordering must run fastest in contrast to the Quad and Hex ordering where 'p' index runs fastest to be consistent with the quadrature ordering.

In the triangular space the i (i.e. $\eta_1$ direction) ordering still runs fastest by convention.

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::NodalTriExp, Nektar::LocalRegions::TriExp, and Nektar::StdRegions::StdNodalTriExp.

Definition at line 460 of file StdTriExp.cpp.

References v_IProductWRTBase_SumFac().

Referenced by v_FwdTrans(), and v_FwdTrans_BndConstrained().

463  {
464  StdTriExp::v_IProductWRTBase_SumFac(inarray,outarray);
465  }
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
Definition: StdTriExp.cpp:479
void Nektar::StdRegions::StdTriExp::v_IProductWRTBase_MatOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented in Nektar::LocalRegions::TriExp.

Definition at line 467 of file StdTriExp.cpp.

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

470  {
471  int nq = GetTotPoints();
472  StdMatrixKey iprodmatkey(eIProductWRTBase,DetShapeType(),*this);
473  DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
474 
475  Blas::Dgemv('N',m_ncoeffs,nq,1.0,iprodmat->GetPtr().get(),
476  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
477  }
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:706
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:141
void Nektar::StdRegions::StdTriExp::v_IProductWRTBase_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
bool  multiplybyweights = true 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::NodalTriExp, Nektar::LocalRegions::TriExp, and Nektar::StdRegions::StdNodalTriExp.

Definition at line 479 of file StdTriExp.cpp.

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

Referenced by v_IProductWRTBase(), and Nektar::StdRegions::StdNodalTriExp::v_IProductWRTBase_SumFac().

483  {
484  int nquad0 = m_base[0]->GetNumPoints();
485  int nquad1 = m_base[1]->GetNumPoints();
486  int order0 = m_base[0]->GetNumModes();
487 
488  if(multiplybyweights)
489  {
490  Array<OneD,NekDouble> tmp(nquad0*nquad1+nquad1*order0);
491  Array<OneD,NekDouble> wsp(tmp+nquad0*nquad1);
492 
493  // multiply by integration constants
494  MultiplyByQuadratureMetric(inarray,tmp);
495  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
496  m_base[1]->GetBdata(),
497  tmp,outarray,wsp);
498  }
499  else
500  {
501  Array<OneD,NekDouble> wsp(nquad1*order0);
502  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
503  m_base[1]->GetBdata(),
504  inarray,outarray,wsp);
505  }
506  }
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:947
void IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdTriExp::v_IProductWRTBase_SumFacKernel ( const Array< OneD, const NekDouble > &  base0,
const Array< OneD, const NekDouble > &  base1,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
Array< OneD, NekDouble > &  wsp,
bool  doCheckCollDir0,
bool  doCheckCollDir1 
)
protectedvirtual

Implements Nektar::StdRegions::StdExpansion2D.

Definition at line 508 of file StdTriExp.cpp.

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

516  {
517  int i;
518  int mode;
519  int nquad0 = m_base[0]->GetNumPoints();
520  int nquad1 = m_base[1]->GetNumPoints();
521  int nmodes0 = m_base[0]->GetNumModes();
522  int nmodes1 = m_base[1]->GetNumModes();
523 
524  ASSERTL1(wsp.num_elements() >= nquad1*nmodes0,
525  "Workspace size is not sufficient");
526 
527  Blas::Dgemm('T','N',nquad1,nmodes0,nquad0,1.0,inarray.get(),nquad0,
528  base0.get(),nquad0,0.0,wsp.get(),nquad1);
529 
530  // Inner product with respect to 'b' direction
531  for (mode=i=0; i < nmodes0; ++i)
532  {
533  Blas::Dgemv('T',nquad1,nmodes1-i,1.0, base1.get()+mode*nquad1,
534  nquad1,wsp.get()+i*nquad1,1, 0.0,
535  outarray.get() + mode,1);
536  mode += nmodes1 - i;
537  }
538 
539  // fix for modified basis by splitting top vertex mode
541  {
542  outarray[1] += Blas::Ddot(nquad1,base1.get()+nquad1,1,
543  wsp.get()+nquad1,1);
544  }
545  }
Principle Modified Functions .
Definition: BasisType.h:49
T Ddot(int n, const Array< OneD, const T > &w, const int incw, const Array< OneD, const T > &x, const int incx, const Array< OneD, const int > &y, const int incy)
Definition: VmathArray.hpp:436
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:228
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdTriExp::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::NodalTriExp, Nektar::StdRegions::StdNodalTriExp, and Nektar::LocalRegions::TriExp.

Definition at line 547 of file StdTriExp.cpp.

References v_IProductWRTDerivBase_SumFac().

551  {
552  StdTriExp::v_IProductWRTDerivBase_SumFac(dir,inarray,outarray);
553  }
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdTriExp.cpp:589
void Nektar::StdRegions::StdTriExp::v_IProductWRTDerivBase_MatOp ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented in Nektar::LocalRegions::TriExp.

Definition at line 555 of file StdTriExp.cpp.

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

559  {
560  int nq = GetTotPoints();
562 
563  switch(dir)
564  {
565  case 0:
566  {
567  mtype = eIProductWRTDerivBase0;
568  break;
569  }
570  case 1:
571  {
572  mtype = eIProductWRTDerivBase1;
573  break;
574  }
575  default:
576  {
577  ASSERTL1(false,"input dir is out of range");
578  break;
579  }
580  }
581 
582  StdMatrixKey iprodmatkey(mtype,DetShapeType(),*this);
583  DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
584 
585  Blas::Dgemv('N',m_ncoeffs,nq,1.0,iprodmat->GetPtr().get(),
586  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
587  }
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:706
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:141
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
void Nektar::StdRegions::StdTriExp::v_IProductWRTDerivBase_SumFac ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::NodalTriExp, Nektar::LocalRegions::TriExp, and Nektar::StdRegions::StdNodalTriExp.

Definition at line 589 of file StdTriExp.cpp.

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

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

593  {
594  int i;
595  int nquad0 = m_base[0]->GetNumPoints();
596  int nquad1 = m_base[1]->GetNumPoints();
597  int nqtot = nquad0*nquad1;
598  int nmodes0 = m_base[0]->GetNumModes();
599  int wspsize = max(max(nqtot,m_ncoeffs),nquad1*nmodes0);
600 
601  Array<OneD, NekDouble> gfac0(2*wspsize);
602  Array<OneD, NekDouble> tmp0 (gfac0+wspsize);
603 
604  const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
605 
606  // set up geometric factor: 2/(1-z1)
607  for(i = 0; i < nquad1; ++i)
608  {
609  gfac0[i] = 2.0/(1-z1[i]);
610  }
611 
612  for(i = 0; i < nquad1; ++i)
613  {
614  Vmath::Smul(nquad0,gfac0[i],&inarray[0]+i*nquad0,1,
615  &tmp0[0]+i*nquad0,1);
616  }
617 
618  MultiplyByQuadratureMetric(tmp0,tmp0);
619 
620  switch(dir)
621  {
622  case 0:
623  {
624  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
625  m_base[1]->GetBdata(),
626  tmp0,outarray,gfac0);
627  break;
628  }
629  case 1:
630  {
631  Array<OneD, NekDouble> tmp3(m_ncoeffs);
632  const Array<OneD, const NekDouble>& z0 = m_base[0]->GetZ();
633 
634  for (i = 0; i < nquad0; ++i)
635  {
636  gfac0[i] = 0.5*(1+z0[i]);
637  }
638 
639  for (i = 0; i < nquad1; ++i)
640  {
641  Vmath::Vmul(nquad0,&gfac0[0],1,&tmp0[0]+i*nquad0,1,
642  &tmp0[0]+i*nquad0,1);
643  }
644 
645  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
646  m_base[1]->GetBdata(),
647  tmp0,tmp3,gfac0);
648 
649  MultiplyByQuadratureMetric(inarray,tmp0);
650  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
651  m_base[1]->GetDbdata(),
652  tmp0,outarray,gfac0);
653  Vmath::Vadd(m_ncoeffs,&tmp3[0],1,&outarray[0],1,
654  &outarray[0],1);
655  break;
656  }
657  default:
658  {
659  ASSERTL1(false, "input dir is out of range");
660  break;
661  }
662  }
663  }
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:947
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:213
void IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:299
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:183
bool Nektar::StdRegions::StdTriExp::v_IsBoundaryInteriorExpansion ( )
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 845 of file StdTriExp.cpp.

References Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, and Nektar::StdRegions::StdExpansion::m_base.

846  {
847  return m_base[0]->GetBasisType() == LibUtilities::eModified_A &&
848  m_base[1]->GetBasisType() == LibUtilities::eModified_B;
849  }
Principle Modified Functions .
Definition: BasisType.h:49
Principle Modified Functions .
Definition: BasisType.h:50
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdTriExp::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::NodalTriExp, Nektar::LocalRegions::TriExp, and Nektar::StdRegions::StdNodalTriExp.

Definition at line 1347 of file StdTriExp.cpp.

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

1351  {
1352  StdTriExp::v_LaplacianMatrixOp_MatFree(inarray,outarray,mkey);
1353  }
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
void Nektar::StdRegions::StdTriExp::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::NodalTriExp, Nektar::LocalRegions::TriExp, and Nektar::StdRegions::StdNodalTriExp.

Definition at line 1355 of file StdTriExp.cpp.

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

1361  {
1363  k1,k2,inarray,outarray,mkey);
1364  }
void LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void Nektar::StdRegions::StdTriExp::v_LocCoordToLocCollapsed ( const Array< OneD, const NekDouble > &  xi,
Array< OneD, NekDouble > &  eta 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 669 of file StdTriExp.cpp.

References Nektar::NekConstants::kNekZeroTol.

671  {
672 
673  // set up local coordinate system
674  if (fabs(xi[1]-1.0) < NekConstants::kNekZeroTol)
675  {
676  eta[0] = -1.0;
677  eta[1] = 1.0;
678  }
679  else
680  {
681  eta[0] = 2*(1+xi[0])/(1-xi[1])-1.0;
682  eta[1] = xi[1];
683  }
684  }
static const NekDouble kNekZeroTol
void Nektar::StdRegions::StdTriExp::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::NodalTriExp, Nektar::LocalRegions::TriExp, and Nektar::StdRegions::StdNodalTriExp.

Definition at line 1339 of file StdTriExp.cpp.

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

1343  {
1344  StdExpansion::MassMatrixOp_MatFree(inarray,outarray,mkey);
1345  }
void MassMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void Nektar::StdRegions::StdTriExp::v_MultiplyByStdQuadratureMetric ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1559 of file StdTriExp.cpp.

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

1562  {
1563  int i;
1564  int nquad0 = m_base[0]->GetNumPoints();
1565  int nquad1 = m_base[1]->GetNumPoints();
1566 
1567  const Array<OneD, const NekDouble>& w0 = m_base[0]->GetW();
1568  const Array<OneD, const NekDouble>& w1 = m_base[1]->GetW();
1569  const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
1570 
1571  // multiply by integration constants
1572  for(i = 0; i < nquad1; ++i)
1573  {
1574  Vmath::Vmul(nquad0,inarray.get()+i*nquad0,1,
1575  w0.get(),1, outarray.get()+i*nquad0,1);
1576  }
1577 
1578  switch(m_base[1]->GetPointsType())
1579  {
1580  // Legendre inner product
1583  for(i = 0; i < nquad1; ++i)
1584  {
1585  Blas::Dscal(nquad0,0.5*(1-z1[i])*w1[i],
1586  outarray.get()+i*nquad0,1);
1587  }
1588  break;
1589 
1590  // (1,0) Jacobi Inner product
1592  for(i = 0; i < nquad1; ++i)
1593  {
1594  Blas::Dscal(nquad0,0.5*w1[i],outarray.get()+i*nquad0,1);
1595  }
1596  break;
1597 
1598  default:
1599  ASSERTL0(false, "Unsupported quadrature points type.");
1600  break;
1601  }
1602  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:59
1D Evenly-spaced points using Lagrange polynomial
Definition: PointsType.h:65
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
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:52
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:183
int Nektar::StdRegions::StdTriExp::v_NumBndryCoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp.

Definition at line 750 of file StdTriExp.cpp.

References ASSERTL1, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::StdRegions::StdExpansion::GetBasisNumModes(), and Nektar::StdRegions::StdExpansion::GetBasisType().

Referenced by v_FwdTrans_BndConstrained().

751  {
753  "BasisType is not a boundary interior form");
755  "BasisType is not a boundary interior form");
756 
757  return 3 + (GetBasisNumModes(0)-2) + 2*(GetBasisNumModes(1)-2);
758  }
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:178
Principle Modified Functions .
Definition: BasisType.h:49
Principle Modified Functions .
Definition: BasisType.h:50
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:228
int Nektar::StdRegions::StdTriExp::v_NumDGBndryCoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 760 of file StdTriExp.cpp.

References ASSERTL1, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::StdRegions::StdExpansion::GetBasisNumModes(), and Nektar::StdRegions::StdExpansion::GetBasisType().

761  {
763  "BasisType is not a boundary interior form");
765  "BasisType is not a boundary interior form");
766 
767  return GetBasisNumModes(0) + 2*GetBasisNumModes(1);
768  }
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:178
Principle Modified Functions .
Definition: BasisType.h:49
Principle Modified Functions .
Definition: BasisType.h:50
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:228
void Nektar::StdRegions::StdTriExp::v_PhysDeriv ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d0,
Array< OneD, NekDouble > &  out_d1,
Array< OneD, NekDouble > &  out_d2 = NullNekDouble1DArray 
)
protectedvirtual

Calculate the derivative of the physical points.

$ \frac{\partial u}{\partial x_1} = \left . \frac{2.0}{1-\eta_2} \frac{\partial u}{\partial d\eta_1} \right |_{\eta_2}$

$ \frac{\partial u}{\partial x_2} = \left . \frac{1+\eta_1}{1-\eta_2} \frac{\partial u}{\partial d\eta_1} \right |_{\eta_2} + \left . \frac{\partial u}{\partial d\eta_2} \right |_{\eta_1} $

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::NodalTriExp, and Nektar::LocalRegions::TriExp.

Definition at line 130 of file StdTriExp.cpp.

References Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion2D::PhysTensorDeriv(), Vmath::Sadd(), Vmath::Sdiv(), Vmath::Smul(), and Vmath::Vvtvp().

Referenced by v_PhysDeriv(), and v_StdPhysDeriv().

135  {
136  int i;
137  int nquad0 = m_base[0]->GetNumPoints();
138  int nquad1 = m_base[1]->GetNumPoints();
139  Array<OneD, NekDouble> wsp(std::max(nquad0, nquad1));
140 
141  const Array<OneD, const NekDouble>& z0 = m_base[0]->GetZ();
142  const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
143 
144  // set up geometric factor: 2/(1-z1)
145  Vmath::Sadd(nquad1, -1.0, z1, 1, wsp, 1);
146  Vmath::Sdiv(nquad1, -2.0, wsp, 1, wsp, 1);
147 
148  if (out_d0.num_elements() > 0)
149  {
150  PhysTensorDeriv(inarray, out_d0, out_d1);
151 
152  for (i = 0; i < nquad1; ++i)
153  {
154  Blas::Dscal(nquad0,wsp[i],&out_d0[0]+i*nquad0,1);
155  }
156 
157  // if no d1 required do not need to calculate both deriv
158  if (out_d1.num_elements() > 0)
159  {
160  // set up geometric factor: (1_z0)/(1-z1)
161  Vmath::Sadd(nquad0, 1.0, z0, 1, wsp, 1);
162  Vmath::Smul(nquad0, 0.5, wsp, 1, wsp, 1);
163 
164  for (i = 0; i < nquad1; ++i)
165  {
166  Vmath::Vvtvp(nquad0,&wsp[0],1,&out_d0[0]+i*nquad0,
167  1,&out_d1[0]+i*nquad0,
168  1,&out_d1[0]+i*nquad0,1);
169  }
170  }
171  }
172  else if (out_d1.num_elements() > 0)
173  {
174  Array<OneD, NekDouble> diff0(nquad0*nquad1);
175  PhysTensorDeriv(inarray, diff0, out_d1);
176 
177  for (i = 0; i < nquad1; ++i)
178  {
179  Blas::Dscal(nquad0,wsp[i],&diff0[0]+i*nquad0,1);
180  }
181 
182  Vmath::Sadd(nquad0, 1.0, z0, 1, wsp, 1);
183  Vmath::Smul(nquad0, 0.5, wsp, 1, wsp, 1);
184 
185  for (i = 0; i < nquad1; ++i)
186  {
187  Vmath::Vvtvp(nquad0,&wsp[0],1,&diff0[0]+i*nquad0,
188  1,&out_d1[0]+i*nquad0,
189  1,&out_d1[0]+i*nquad0,1);
190  }
191  }
192  }
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:442
void Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/y.
Definition: Vmath.cpp:271
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d0, Array< OneD, NekDouble > &outarray_d1)
Calculate the 2D derivative in the local tensor/collapsed coordinate at the physical points...
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:213
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition: Vmath.cpp:315
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdTriExp::v_PhysDeriv ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d0 
)
protectedvirtual

Calculate the derivative of the physical points in a given direction.

See also
StdRegions::StdExpansion::PhysDeriv

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::NodalTriExp, and Nektar::LocalRegions::TriExp.

Definition at line 194 of file StdTriExp.cpp.

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

198  {
199  switch(dir)
200  {
201  case 0:
202  {
203  v_PhysDeriv(inarray, outarray, NullNekDouble1DArray);
204  break;
205  }
206  case 1:
207  {
208  v_PhysDeriv(inarray, NullNekDouble1DArray, outarray);
209  break;
210  }
211  default:
212  {
213  ASSERTL1(false,"input dir is out of range");
214  break;
215  }
216  }
217  }
static Array< OneD, NekDouble > NullNekDouble1DArray
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=NullNekDouble1DArray)
Calculate the derivative of the physical points.
Definition: StdTriExp.cpp:130
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
void Nektar::StdRegions::StdTriExp::v_ReduceOrderCoeffs ( int  numMin,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp.

Definition at line 1476 of file StdTriExp.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::eOrtho_B, Nektar::StdRegions::StdExpansion::GetBasisType(), and Nektar::StdRegions::StdExpansion::m_base.

1480  {
1481  int n_coeffs = inarray.num_elements();
1482  int nquad0 = m_base[0]->GetNumPoints();
1483  int nquad1 = m_base[1]->GetNumPoints();
1484  Array<OneD, NekDouble> coeff(n_coeffs);
1485  Array<OneD, NekDouble> coeff_tmp(n_coeffs,0.0);
1486  Array<OneD, NekDouble> tmp;
1487  Array<OneD, NekDouble> tmp2;
1488  int nqtot = nquad0*nquad1;
1489  Array<OneD, NekDouble> phys_tmp(nqtot,0.0);
1490 
1491  int nmodes0 = m_base[0]->GetNumModes();
1492  int nmodes1 = m_base[1]->GetNumModes();
1493  int numMin2 = nmodes0;
1494  int i;
1495 
1496  const LibUtilities::PointsKey Pkey0(
1498  const LibUtilities::PointsKey Pkey1(
1500 
1501  LibUtilities::BasisKey b0(m_base[0]->GetBasisType(),nmodes0,Pkey0);
1502  LibUtilities::BasisKey b1(m_base[1]->GetBasisType(),nmodes1,Pkey1);
1503 
1504  LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A,nmodes0,Pkey0);
1505  LibUtilities::BasisKey bortho1(LibUtilities::eOrtho_B,nmodes1,Pkey1);
1506 
1507  StdRegions::StdTriExpSharedPtr m_OrthoTriExp;
1509 
1511  ::AllocateSharedPtr(b0, b1);
1512  m_OrthoTriExp = MemoryManager<StdRegions::StdTriExp>
1513  ::AllocateSharedPtr(bortho0, bortho1);
1514 
1515  m_TriExp ->BwdTrans(inarray,phys_tmp);
1516  m_OrthoTriExp->FwdTrans(phys_tmp, coeff);
1517 
1518  for (i = 0; i < n_coeffs; i++)
1519  {
1520  if (i == numMin)
1521  {
1522  coeff[i] = 0.0;
1523  numMin += numMin2 - 1;
1524  numMin2 -= 1.0;
1525  }
1526  }
1527 
1528  m_OrthoTriExp->BwdTrans(coeff,phys_tmp);
1529  m_TriExp ->FwdTrans(phys_tmp, outarray);
1530 
1531  }
boost::shared_ptr< StdTriExp > StdTriExpSharedPtr
Definition: StdTriExp.h:267
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
Principle Orthogonal Functions .
Definition: BasisType.h:47
Principle Orthogonal Functions .
Definition: BasisType.h:46
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
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:52
void Nektar::StdRegions::StdTriExp::v_StdPhysDeriv ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d0,
Array< OneD, NekDouble > &  out_d1,
Array< OneD, NekDouble > &  out_d2 = NullNekDouble1DArray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::NodalTriExp.

Definition at line 219 of file StdTriExp.cpp.

References v_PhysDeriv().

224  {
225  StdTriExp::v_PhysDeriv(inarray, out_d0, out_d1);
226  }
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=NullNekDouble1DArray)
Calculate the derivative of the physical points.
Definition: StdTriExp.cpp:130
void Nektar::StdRegions::StdTriExp::v_StdPhysDeriv ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 228 of file StdTriExp.cpp.

References v_PhysDeriv().

232  {
233  StdTriExp::v_PhysDeriv(dir,inarray,outarray);
234  }
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=NullNekDouble1DArray)
Calculate the derivative of the physical points.
Definition: StdTriExp.cpp:130
void Nektar::StdRegions::StdTriExp::v_SVVLaplacianFilter ( Array< OneD, NekDouble > &  array,
const StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp.

Definition at line 1384 of file StdTriExp.cpp.

References Nektar::StdRegions::StdExpansion::BwdTrans(), Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::eFactorSVVDiffCoeff, Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::eOrtho_B, Nektar::StdRegions::eVarCoeffLaplacian, Nektar::StdRegions::StdExpansion::FwdTrans(), Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::StdRegions::StdMatrixKey::GetVarCoeff(), Nektar::StdRegions::StdMatrixKey::HasVarCoeff(), Nektar::StdRegions::StdExpansion::m_base, Vmath::Vmul(), and Vmath::Vsqrt().

1386  {
1387  int qa = m_base[0]->GetNumPoints();
1388  int qb = m_base[1]->GetNumPoints();
1389  int nmodes_a = m_base[0]->GetNumModes();
1390  int nmodes_b = m_base[1]->GetNumModes();
1391  int nmodes = min(nmodes_a,nmodes_b);
1392 
1393  // Declare orthogonal basis.
1394  LibUtilities::PointsKey pa(qa,m_base[0]->GetPointsType());
1395  LibUtilities::PointsKey pb(qb,m_base[1]->GetPointsType());
1396 
1397  LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A,nmodes_a,pa);
1398  LibUtilities::BasisKey Bb(LibUtilities::eOrtho_B,nmodes_b,pb);
1399  StdTriExp OrthoExp(Ba,Bb);
1400 
1401  NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDiffCoeff);
1402 
1403  Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
1404 
1405 
1406  if(mkey.HasVarCoeff(eVarCoeffLaplacian)) // Rodrigo's svv mapping
1407  {
1408  Array<OneD, NekDouble> sqrt_varcoeff(qa*qb);
1409  Array<OneD, NekDouble> tmp(qa*qb);
1410 
1411  Vmath::Vsqrt(qa * qb,
1412  mkey.GetVarCoeff(eVarCoeffLaplacian), 1,
1413  sqrt_varcoeff, 1);
1414 
1415  // multiply by sqrt(Variable Coefficient) containing h v /p
1416  Vmath::Vmul(qa*qb,sqrt_varcoeff,1,array,1,tmp,1);
1417 
1418  // project onto modal space.
1419  OrthoExp.FwdTrans(tmp,orthocoeffs);
1420 
1421  int cnt = 0;
1422  for(int j = 0; j < nmodes_a; ++j)
1423  {
1424  for(int k = 0; k < nmodes_b-j; ++k, ++cnt)
1425  {
1426  orthocoeffs[cnt] *=
1427  (1.0 + SvvDiffCoeff
1428  *pow(j/(nmodes_a-1)+k/(nmodes_b-1),0.5*nmodes));
1429  }
1430  }
1431 
1432  // backward transform to physical space
1433  OrthoExp.BwdTrans(orthocoeffs,tmp);
1434 
1435  // multiply by sqrt(Variable Coefficient) containing h v /p
1436  // - split to keep symmetry
1437  Vmath::Vmul(qa*qb,sqrt_varcoeff,1,tmp,1,array,1);
1438  }
1439  else
1440  {
1441  int j, k , cnt = 0;
1442  int cutoff = (int) (mkey.GetConstFactor(eFactorSVVCutoffRatio)*
1443  min(nmodes_a,nmodes_b));
1444 
1445  NekDouble epsilon = 1.0;
1446  int nmodes = min(nmodes_a,nmodes_b);
1447 
1448  // project onto physical space.
1449  OrthoExp.FwdTrans(array,orthocoeffs);
1450 
1451  // apply SVV filter (JEL)
1452  for(j = 0; j < nmodes_a; ++j)
1453  {
1454  for(k = 0; k < nmodes_b-j; ++k)
1455  {
1456  if(j + k >= cutoff)
1457  {
1458  orthocoeffs[cnt] *= (SvvDiffCoeff
1459  *exp(-(j+k-nmodes)*(j+k-nmodes)
1460  /((NekDouble)((j+k-cutoff+epsilon)
1461  *(j+k-cutoff+epsilon)))));
1462  }
1463  else
1464  {
1465  orthocoeffs[cnt] *= 0.0;
1466  }
1467  cnt++;
1468  }
1469  }
1470 
1471  // backward transform to physical space
1472  OrthoExp.BwdTrans(orthocoeffs,array);
1473  }
1474  }
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
Definition: Vmath.cpp:408
Principle Orthogonal Functions .
Definition: BasisType.h:47
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 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:183
void Nektar::StdRegions::StdTriExp::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::NodalTriExp, Nektar::LocalRegions::TriExp, and Nektar::StdRegions::StdNodalTriExp.

Definition at line 1366 of file StdTriExp.cpp.

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

1371  {
1372  StdExpansion::WeakDerivMatrixOp_MatFree(i,inarray,outarray,mkey);
1373  }
void WeakDerivMatrixOp_MatFree(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)