Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Protected Member Functions | List of all members
Nektar::StdRegions::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
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion ()
 
bool IsNodalNonTensorialExp ()
 
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space. More...
 
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space. More...
 
void FwdTrans_BndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble Integral (const Array< OneD, const NekDouble > &inarray)
 This function integrates the specified function over the domain. More...
 
void FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 This function fills the array outarray with the mode-th mode of the expansion. More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 this function calculates the inner product of a given function f with the different modes of the expansion More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check)
 
void IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
int GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void SetElmtId (const int id)
 Set the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2=NullNekDouble1DArray, Array< OneD, NekDouble > &coords_3=NullNekDouble1DArray)
 this function returns the physical coordinates of the quadrature points of the expansion More...
 
void GetCoord (const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
 given the coordinates of a point of the element in the local collapsed coordinate system, this function calculates the physical coordinates of the point More...
 
DNekMatSharedPtr GetStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr GetStdStaticCondMatrix (const StdMatrixKey &mkey)
 
IndexMapValuesSharedPtr GetIndexMap (const IndexMapKey &ikey)
 
const Array< OneD, const
NekDouble > & 
GetPhysNormals (void)
 
void SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
virtual void SetUpPhysNormals (const int edge)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
void DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
StdRegions::Orientation GetForient (int face)
 
StdRegions::Orientation GetEorient (int edge)
 
StdRegions::Orientation GetPorient (int point)
 
StdRegions::Orientation GetCartesianEorient (int edge)
 
void SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
 
void SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
int CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
void ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs)
 
NekDouble StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
int GetCoordim ()
 
void GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
void GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
int GetVertexMap (const int localVertexId, bool useCoeffPacking=false)
 
void GetEdgeInteriorMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
 
void GetFaceInteriorMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
 
void GetEdgeToElementMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int P=-1)
 
void GetFaceToElementMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int nummodesA=-1, int nummodesB=-1)
 
void GetEdgePhysVals (const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Extract the physical values along edge edge from inarray into outarray following the local edge orientation and point distribution defined by defined in EdgeExp. More...
 
void GetEdgePhysVals (const int edge, const boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetTracePhysVals (const int edge, const boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetVertexPhysVals (const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray)
 
void GetEdgeInterpVals (const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetEdgeQFactors (const int edge, Array< OneD, NekDouble > &outarray)
 Extract the metric factors to compute the contravariant fluxes along edge edge and stores them into outarray following the local edge orientation (i.e. anticlockwise convention). More...
 
void GetFacePhysVals (const int face, const boost::shared_ptr< StdExpansion > &FaceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=eNoOrientation)
 
void GetEdgePhysMap (const int edge, Array< OneD, int > &outarray)
 
void GetFacePhysMap (const int face, Array< OneD, int > &outarray)
 
void MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
DNekMatSharedPtr CreateGeneralMatrix (const StdMatrixKey &mkey)
 this function generates the mass matrix $\mathbf{M}[i][j] = \int \phi_i(\mathbf{x}) \phi_j(\mathbf{x}) d\mathbf{x}$ More...
 
void GeneralMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
DNekMatSharedPtr GenMatrix (const StdMatrixKey &mkey)
 
void PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void PhysDeriv_s (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds)
 
void PhysDeriv_n (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn)
 
void PhysDirectionalDeriv (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &outarray)
 
void StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void AddRobinMassMatrix (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
void AddRobinEdgeContribution (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, Array< OneD, NekDouble > &coeffs)
 
NekDouble PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta. More...
 
const boost::shared_ptr
< SpatialDomains::GeomFactors > & 
GetMetricInfo (void) const
 
virtual int v_GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
virtual const Array< OneD,
const NekDouble > & 
v_GetPhysNormals (void)
 
virtual void v_SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
virtual void v_SetUpPhysNormals (const int edge)
 
virtual void v_ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmode_offset, NekDouble *coeffs)
 Unpack data from input file assuming it comes from the same expansion type. More...
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual void v_DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual StdRegions::Orientation v_GetForient (int face)
 
virtual StdRegions::Orientation v_GetEorient (int edge)
 
virtual StdRegions::Orientation v_GetCartesianEorient (int edge)
 
virtual StdRegions::Orientation v_GetPorient (int point)
 
NekDouble Linf (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete $ L_\infty$ error $ |\epsilon|_\infty = \max |u - u_{exact}|$ where $ u_{exact}$ is given by the array sol. More...
 
NekDouble L2 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete $ L_2$ error, $ | \epsilon |_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 dx \right]^{1/2} d\xi_1 $ where $ u_{exact}$ is given by the array sol. More...
 
NekDouble H1 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete $ H^1$ error, $ | \epsilon |^1_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 + \nabla(u - u_{exact})\cdot\nabla(u - u_{exact})\cdot dx \right]^{1/2} d\xi_1 $ where $ u_{exact}$ is given by the array sol. More...
 
const NormalVectorGetEdgeNormal (const int edge) const
 
void ComputeEdgeNormal (const int edge)
 
void NegateEdgeNormal (const int edge)
 
bool EdgeNormalNegated (const int edge)
 
void ComputeFaceNormal (const int face)
 
void NegateFaceNormal (const int face)
 
bool FaceNormalNegated (const int face)
 
void ComputeVertexNormal (const int vertex)
 
const NormalVectorGetFaceNormal (const int face) const
 
const NormalVectorGetVertexNormal (const int vertex) const
 
const NormalVectorGetSurfaceNormal (const int id) const
 
const LibUtilities::PointsKeyVector GetPointsKeys () const
 
Array< OneD, unsigned int > GetEdgeInverseBoundaryMap (int eid)
 
Array< OneD, unsigned int > GetFaceInverseBoundaryMap (int fid, StdRegions::Orientation faceOrient=eNoOrientation)
 
DNekMatSharedPtr BuildInverseTransformationMatrix (const DNekScalMatSharedPtr &m_transformationmatrix)
 
void PhysInterpToSimplexEquiSpaced (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int npset=-1)
 This function performs an interpolation from the physical space points provided at input into an array of equispaced points which are not the collapsed coordinate. So for a tetrahedron you will only get a tetrahedral number of values. More...
 
void GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true)
 This function provides the connectivity of local simplices (triangles or tets) to connect the equispaced data points provided by PhysInterpToSimplexEquiSpaced. More...
 
void EquiSpacedToCoeffs (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs a projection/interpolation from the equispaced points sometimes used in post-processing onto the coefficient space. More...
 
template<class T >
boost::shared_ptr< T > as ()
 
void IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 

Protected Member Functions

virtual 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 45 of file StdTriExp.cpp.

46  {
47  }
Nektar::StdRegions::StdTriExp::StdTriExp ( const LibUtilities::BasisKey Ba,
const LibUtilities::BasisKey Bb 
)

Definition at line 50 of file StdTriExp.cpp.

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

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

Definition at line 67 of file StdTriExp.cpp.

67  :
68  StdExpansion(T),
70  {
71  }
StdExpansion()
Default Constructor.
Nektar::StdRegions::StdTriExp::~StdTriExp ( )

Definition at line 73 of file StdTriExp.cpp.

74  {
75  }

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 250 of file StdTriExp.cpp.

References v_BwdTrans_SumFac().

253  {
254  v_BwdTrans_SumFac(inarray,outarray);
255  }
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdTriExp.cpp:258
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 258 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().

261  {
262  Array<OneD, NekDouble> wsp(m_base[0]->GetNumPoints()*
263  m_base[1]->GetNumModes());
264 
265  BwdTrans_SumFacKernel(m_base[0]->GetBdata(),
266  m_base[1]->GetBdata(),
267  inarray,outarray,wsp);
268  }
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 270 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.

278  {
279  int i;
280  int mode;
281  int nquad0 = m_base[0]->GetNumPoints();
282  int nquad1 = m_base[1]->GetNumPoints();
283  int nmodes0 = m_base[0]->GetNumModes();
284  int nmodes1 = m_base[1]->GetNumModes();
285 
286  ASSERTL1(wsp.num_elements() >= nquad0*nmodes1,
287  "Workspace size is not sufficient");
290  "Basis[1] is not of general tensor type");
291 
292  for (i = mode = 0; i < nmodes0; ++i)
293  {
294  Blas::Dgemv('N', nquad1,nmodes1-i,1.0,base1.get()+mode*nquad1,
295  nquad1,&inarray[0]+mode,1,0.0,&wsp[0]+i*nquad1,1);
296  mode += nmodes1-i;
297  }
298 
299  // fix for modified basis by splitting top vertex mode
301  {
302  Blas::Daxpy(nquad1,inarray[1],base1.get()+nquad1,1,
303  &wsp[0]+nquad1,1);
304  }
305 
306  Blas::Dgemm('N','T', nquad0,nquad1,nmodes0,1.0, base0.get(),nquad0,
307  &wsp[0], nquad1,0.0, &outarray[0], nquad0);
308  }
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:213
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
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 802 of file StdTriExp.cpp.

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

805  {
807  nummodes[modes_offset],
808  nummodes[modes_offset+1]);
809  modes_offset += 2;
810 
811  return nmodes;
812  }
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 1333 of file StdTriExp.cpp.

References v_GenMatrix().

1334  {
1335  return v_GenMatrix(mkey);
1336  }
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey)
Definition: StdTriExp.cpp:1256
int Nektar::StdRegions::StdTriExp::v_DetCartesianDirOfEdge ( const int  edge)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 855 of file StdTriExp.cpp.

References ASSERTL2.

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 862 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().

864  {
865  ASSERTL2(i >= 0 && i <= 2, "edge id is out of range");
866 
867  if (i == 0)
868  {
869  return GetBasis(0)->GetBasisKey();
870  }
871  else
872  {
873  switch(m_base[1]->GetBasisType())
874  {
876  {
877  switch(m_base[1]->GetPointsType())
878  {
880  {
881  LibUtilities::PointsKey pkey(
882  m_base[1]->GetBasisKey().GetPointsKey().
883  GetNumPoints()+1,
885  return LibUtilities::BasisKey(
887  m_base[1]->GetNumModes(),pkey);
888  break;
889  }
890 
891  default:
892  ASSERTL0(false,"unexpected points distribution");
893  break;
894  }
895  }
896  default:
897  ASSERTL0(false,"Information not available to set edge key");
898  break;
899  }
900  }
902  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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:57
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:213
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:50
LibUtilities::ShapeType Nektar::StdRegions::StdTriExp::v_DetShapeType ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 749 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 690 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().

692  {
693  int i,m;
694  int nquad0 = m_base[0]->GetNumPoints();
695  int nquad1 = m_base[1]->GetNumPoints();
696  int order0 = m_base[0]->GetNumModes();
697  int order1 = m_base[1]->GetNumModes();
698  int mode0 = 0;
699  Array<OneD, const NekDouble> base0 = m_base[0]->GetBdata();
700  Array<OneD, const NekDouble> base1 = m_base[1]->GetBdata();
701 
702  ASSERTL2(mode <= m_ncoeffs,
703  "calling argument mode is larger than "
704  "total expansion order");
705 
706  m = order1;
707  for (i = 0; i < order0; ++i, m+=order1-i)
708  {
709  if (m > mode)
710  {
711  mode0 = i;
712  break;
713  }
714  }
715 
716  // deal with top vertex mode in modified basis
717  if (mode == 1 &&
719  {
720  Vmath::Fill(nquad0*nquad1 , 1.0, outarray, 1);
721  }
722  else
723  {
724  for (i = 0; i < nquad1; ++i)
725  {
726  Vmath::Vcopy(nquad0,(NekDouble *)(base0.get()+mode0*nquad0),
727  1,&outarray[0]+i*nquad0,1);
728  }
729  }
730 
731  for(i = 0; i < nquad0; ++i)
732  {
733  Vmath::Vmul(nquad1,(NekDouble *)(base1.get() + mode*nquad1),
734  1,&outarray[0]+i,nquad0,&outarray[0]+i,nquad0);
735  }
736  }
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:213
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:169
void Nektar::StdRegions::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 310 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().

313  {
314  v_IProductWRTBase(inarray,outarray);
315 
316  // get Mass matrix inverse
317  StdMatrixKey masskey(eInvMass,DetShapeType(),*this);
318  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
319 
320  // copy inarray in case inarray == outarray
321  NekVector<NekDouble> in(m_ncoeffs,outarray,eCopy);
322  NekVector<NekDouble> out(m_ncoeffs,outarray,eWrapper);
323 
324  out = (*matsys)*in;
325  }
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:464
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:700
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 328 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().

331  {
332  int i,j;
333  int npoints[2] = {m_base[0]->GetNumPoints(),
334  m_base[1]->GetNumPoints()};
335  int nmodes[2] = {m_base[0]->GetNumModes(),
336  m_base[1]->GetNumModes()};
337 
338  fill(outarray.get(), outarray.get()+m_ncoeffs, 0.0 );
339 
340  Array<OneD, NekDouble> physEdge[3];
341  Array<OneD, NekDouble> coeffEdge[3];
342  for(i = 0; i < 3; i++)
343  {
344  physEdge[i] = Array<OneD, NekDouble>(npoints[i!=0]);
345  coeffEdge[i] = Array<OneD, NekDouble>(nmodes[i!=0]);
346  }
347 
348  for(i = 0; i < npoints[0]; i++)
349  {
350  physEdge[0][i] = inarray[i];
351  }
352 
353  for(i = 0; i < npoints[1]; i++)
354  {
355  physEdge[1][i] = inarray[npoints[0]-1+i*npoints[0]];
356  physEdge[2][i] = inarray[(npoints[1]-1)*npoints[0]-i*npoints[0]];
357  }
358 
359  StdSegExpSharedPtr segexp[2] = {
361  m_base[0]->GetBasisKey()),
363  m_base[1]->GetBasisKey())
364  };
365 
366  Array<OneD, unsigned int> mapArray;
367  Array<OneD, int> signArray;
368  NekDouble sign;
369 
370  for (i = 0; i < 3; i++)
371  {
372  //segexp[i!=0]->v_FwdTrans_BndConstrained(physEdge[i],coeffEdge[i]);
373  segexp[i!=0]->FwdTrans_BndConstrained(physEdge[i],coeffEdge[i]);
374 
375  v_GetEdgeToElementMap(i,eForwards,mapArray,signArray);
376  for (j = 0; j < nmodes[i != 0]; j++)
377  {
378  sign = (NekDouble) signArray[j];
379  outarray[ mapArray[j] ] = sign * coeffEdge[i][j];
380  }
381  }
382 
383  Array<OneD, NekDouble> tmp0(m_ncoeffs);
384  Array<OneD, NekDouble> tmp1(m_ncoeffs);
385 
386  StdMatrixKey masskey(eMass,DetShapeType(),*this);
387  MassMatrixOp(outarray,tmp0,masskey);
388  v_IProductWRTBase(inarray,tmp1);
389 
390  Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
391 
392  // get Mass matrix inverse (only of interior DOF)
393  // use block (1,1) of the static condensed system
394  // note: this block alreay contains the inverse matrix
395  DNekMatSharedPtr matsys =
396  (m_stdStaticCondMatrixManager[masskey])->GetBlock(1,1);
397 
398  int nBoundaryDofs = v_NumBndryCoeffs();
399  int nInteriorDofs = m_ncoeffs - nBoundaryDofs;
400 
401  Array<OneD, NekDouble> rhs (nInteriorDofs);
402  Array<OneD, NekDouble> result(nInteriorDofs);
403 
404  v_GetInteriorMap(mapArray);
405 
406  for (i = 0; i < nInteriorDofs; i++)
407  {
408  rhs[i] = tmp1[ mapArray[i] ];
409  }
410 
411  Blas::Dgemv('N',nInteriorDofs,nInteriorDofs,
412  1.0,&(matsys->GetPtr())[0],nInteriorDofs,
413  rhs.get(),1,
414  0.0,result.get(),1);
415 
416  for (i = 0; i < nInteriorDofs; i++)
417  {
418  outarray[ mapArray[i] ] = result[i];
419  }
420  }
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:910
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:971
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:22
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:464
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:1189
virtual int v_NumBndryCoeffs() const
Definition: StdTriExp.cpp:754
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:329
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 1537 of file StdTriExp.cpp.

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

1541  {
1543 
1544  if(inarray.get() == outarray.get())
1545  {
1546  Array<OneD,NekDouble> tmp(m_ncoeffs);
1547  Vmath::Vcopy(m_ncoeffs,inarray.get(),1,tmp.get(),1);
1548 
1549  Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, 1.0, mat->GetPtr().get(),
1550  m_ncoeffs, tmp.get(), 1, 0.0, outarray.get(), 1);
1551  }
1552  else
1553  {
1554  Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, 1.0, mat->GetPtr().get(),
1555  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
1556  }
1557  }
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
LibUtilities::NekManager< StdMatrixKey, DNekMat, StdMatrixKey::opLess > m_stdMatrixManager
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
DNekMatSharedPtr Nektar::StdRegions::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 1256 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().

1257  {
1258 
1259  MatrixType mtype = mkey.GetMatrixType();
1260 
1261  DNekMatSharedPtr Mat;
1262 
1263  switch(mtype)
1264  {
1266  {
1267  int nq0, nq1, nq;
1268 
1269  nq0 = m_base[0]->GetNumPoints();
1270  nq1 = m_base[1]->GetNumPoints();
1271 
1272  // take definition from key
1273  if(mkey.ConstFactorExists(eFactorConst))
1274  {
1275  nq = (int) mkey.GetConstFactor(eFactorConst);
1276  }
1277  else
1278  {
1279  nq = max(nq0,nq1);
1280  }
1281 
1282  int neq = LibUtilities::StdTriData::
1284  Array<OneD, Array<OneD, NekDouble> > coords(neq);
1285  Array<OneD, NekDouble> coll (2);
1286  Array<OneD, DNekMatSharedPtr> I (2);
1287  Array<OneD, NekDouble> tmp (nq0);
1288 
1289 
1290  Mat = MemoryManager<DNekMat>::AllocateSharedPtr(neq,nq0*nq1);
1291  int cnt = 0;
1292 
1293  for(int i = 0; i < nq; ++i)
1294  {
1295  for(int j = 0; j < nq-i; ++j,++cnt)
1296  {
1297  coords[cnt] = Array<OneD, NekDouble>(2);
1298  coords[cnt][0] = -1.0 + 2*j/(NekDouble)(nq-1);
1299  coords[cnt][1] = -1.0 + 2*i/(NekDouble)(nq-1);
1300  }
1301  }
1302 
1303  for(int i = 0; i < neq; ++i)
1304  {
1305  LocCoordToLocCollapsed(coords[i],coll);
1306 
1307  I[0] = m_base[0]->GetI(coll);
1308  I[1] = m_base[1]->GetI(coll+1);
1309 
1310  // interpolate first coordinate direction
1311  for (int j = 0; j < nq1; ++j)
1312  {
1313  NekDouble fac = (I[1]->GetPtr())[j];
1314  Vmath::Smul(nq0, fac, I[0]->GetPtr(), 1, tmp, 1);
1315 
1316  Vmath::Vcopy(nq0, &tmp[0], 1,
1317  Mat->GetRawPtr() + j*nq0*neq + i, neq);
1318  }
1319 
1320  }
1321  break;
1322  }
1323  default:
1324  {
1326  break;
1327  }
1328  }
1329 
1330  return Mat;
1331  }
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:199
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:1047
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 1219 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().

1220  {
1223  "Expansion not of expected type");
1224  int i;
1225  int cnt;
1226  int nummodes0, nummodes1;
1227  int value;
1228 
1229  if (outarray.num_elements()!=NumBndryCoeffs())
1230  {
1231  outarray = Array<OneD, unsigned int>(NumBndryCoeffs());
1232  }
1233 
1234  nummodes0 = m_base[0]->GetNumModes();
1235  nummodes1 = m_base[1]->GetNumModes();
1236 
1237  value = 2*nummodes1-1;
1238  for(i = 0; i < value; i++)
1239  {
1240  outarray[i]=i;
1241  }
1242  cnt = value;
1243 
1244  for(i = 0; i < nummodes0-2; i++)
1245  {
1246  outarray[cnt++]=value;
1247  value += nummodes1-2-i;
1248  }
1249  }
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:191
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 829 of file StdTriExp.cpp.

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

832  {
833  Array<OneD, const NekDouble> z0 = m_base[0]->GetZ();
834  Array<OneD, const NekDouble> z1 = m_base[1]->GetZ();
835  int nq0 = GetNumPoints(0);
836  int nq1 = GetNumPoints(1);
837  int i,j;
838 
839  for(i = 0; i < nq1; ++i)
840  {
841  for(j = 0; j < nq0; ++j)
842  {
843  coords_0[i*nq0+j] = (1+z0[j])*(1-z1[i])/2.0 - 1.0;
844  }
845  Vmath::Fill(nq0,z1[i],&coords_1[0] + i*nq0,1);
846  }
847  }
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 814 of file StdTriExp.cpp.

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

815  {
816  ASSERTL2(i >= 0 && i <= 2, "edge id is out of range");
817 
818  if (i == 0)
819  {
820  return GetBasisType(0);
821  }
822  else
823  {
824  return GetBasisType(1);
825  }
826  }
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:213
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 1103 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.

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

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 788 of file StdTriExp.cpp.

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

789  {
790  ASSERTL2((i >= 0)&&(i <= 2),"edge id is out of range");
791 
792  if (i == 0)
793  {
794  return GetNumPoints(0);
795  }
796  else
797  {
798  return GetNumPoints(1);
799  }
800  }
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:213
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 910 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(), and Nektar::StdRegions::StdExpansion::m_base.

Referenced by v_FwdTrans_BndConstrained().

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

1190  {
1193  "Expansion not of a proper type");
1194 
1195  int i,j;
1196  int cnt = 0;
1197  int nummodes0, nummodes1;
1198  int startvalue;
1199  if(outarray.num_elements()!=GetNcoeffs()-NumBndryCoeffs())
1200  {
1201  outarray = Array<OneD, unsigned int>(GetNcoeffs()-NumBndryCoeffs());
1202  }
1203 
1204  nummodes0 = m_base[0]->GetNumModes();
1205  nummodes1 = m_base[1]->GetNumModes();
1206 
1207  startvalue = 2*nummodes1;
1208 
1209  for(i = 0; i < nummodes0-2; i++)
1210  {
1211  for(j = 0; j < nummodes1-3-i; j++)
1212  {
1213  outarray[cnt++]=startvalue+j;
1214  }
1215  startvalue+=nummodes1-2-i;
1216  }
1217  }
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:191
Array< OneD, LibUtilities::BasisSharedPtr > m_base
int Nektar::StdRegions::StdTriExp::v_GetNedges ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 744 of file StdTriExp.cpp.

745  {
746  return 3;
747  }
int Nektar::StdRegions::StdTriExp::v_GetNverts ( ) const
protectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 739 of file StdTriExp.cpp.

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1608 of file StdTriExp.cpp.

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

1611  {
1612  int np1 = m_base[0]->GetNumPoints();
1613  int np2 = m_base[1]->GetNumPoints();
1614  int np = max(np1,np2);
1615 
1616  conn = Array<OneD, int>(3*(np-1)*(np-1));
1617 
1618  int row = 0;
1619  int rowp1 = 0;
1620  int cnt = 0;
1621  for(int i = 0; i < np-1; ++i)
1622  {
1623  rowp1 += np-i;
1624  for(int j = 0; j < np-i-2; ++j)
1625  {
1626  conn[cnt++] = row +j;
1627  conn[cnt++] = row +j+1;
1628  conn[cnt++] = rowp1 +j;
1629 
1630  conn[cnt++] = rowp1 +j+1;
1631  conn[cnt++] = rowp1 +j;
1632  conn[cnt++] = row +j+1;
1633  }
1634 
1635  conn[cnt++] = row +np-i-2;
1636  conn[cnt++] = row +np-i-1;
1637  conn[cnt++] = rowp1+np-i-2;
1638 
1639  row += np-i;
1640  }
1641  }
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 1039 of file StdTriExp.cpp.

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

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

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

1383  {
1384  StdTriExp::v_HelmholtzMatrixOp_MatFree(inarray,outarray,mkey);
1385  }
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 80 of file StdTriExp.cpp.

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

82  {
83  int i;
84  int nquad1 = m_base[1]->GetNumPoints();
85  Array<OneD, NekDouble> w1_tmp(nquad1);
86 
87  Array<OneD, const NekDouble> w0 = m_base[0]->GetW();
88  Array<OneD, const NekDouble> w1 = m_base[1]->GetW();
89  Array<OneD, const NekDouble> z1 = m_base[1]->GetZ();
90 
91  switch(m_base[1]->GetPointsType())
92  {
93  case LibUtilities::eGaussRadauMAlpha1Beta0: // (0,1) Jacobi Inner product
94  {
95  Vmath::Smul(nquad1, 0.5, w1, 1, w1_tmp,1);
96  break;
97  }
98  default:
99  {
100  // include jacobian factor on whatever coordinates are defined.
101  for(i = 0; i < nquad1; ++i)
102  {
103  w1_tmp[i] = 0.5*(1-z1[i])*w1[i];
104  }
105  break;
106  }
107  }
108 
109  return StdExpansion2D::Integral(inarray,w0,w1_tmp);
110  }
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:57
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
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 464 of file StdTriExp.cpp.

References v_IProductWRTBase_SumFac().

Referenced by v_FwdTrans(), and v_FwdTrans_BndConstrained().

467  {
468  StdTriExp::v_IProductWRTBase_SumFac(inarray,outarray);
469  }
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
Definition: StdTriExp.cpp:483
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 471 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.

474  {
475  int nq = GetTotPoints();
476  StdMatrixKey iprodmatkey(eIProductWRTBase,DetShapeType(),*this);
477  DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
478 
479  Blas::Dgemv('N',m_ncoeffs,nq,1.0,iprodmat->GetPtr().get(),
480  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
481  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:470
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:700
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:141
void Nektar::StdRegions::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 483 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().

487  {
488  int nquad0 = m_base[0]->GetNumPoints();
489  int nquad1 = m_base[1]->GetNumPoints();
490  int order0 = m_base[0]->GetNumModes();
491 
492  if(multiplybyweights)
493  {
494  Array<OneD,NekDouble> tmp(nquad0*nquad1+nquad1*order0);
495  Array<OneD,NekDouble> wsp(tmp+nquad0*nquad1);
496 
497  // multiply by integration constants
498  MultiplyByQuadratureMetric(inarray,tmp);
499  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
500  m_base[1]->GetBdata(),
501  tmp,outarray,wsp);
502  }
503  else
504  {
505  Array<OneD,NekDouble> wsp(nquad1*order0);
506  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
507  m_base[1]->GetBdata(),
508  inarray,outarray,wsp);
509  }
510  }
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:942
void IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &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 512 of file StdTriExp.cpp.

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

520  {
521  int i;
522  int mode;
523  int nquad0 = m_base[0]->GetNumPoints();
524  int nquad1 = m_base[1]->GetNumPoints();
525  int nmodes0 = m_base[0]->GetNumModes();
526  int nmodes1 = m_base[1]->GetNumModes();
527 
528  ASSERTL1(wsp.num_elements() >= nquad1*nmodes0,
529  "Workspace size is not sufficient");
530 
531  Blas::Dgemm('T','N',nquad1,nmodes0,nquad0,1.0,inarray.get(),nquad0,
532  base0.get(),nquad0,0.0,wsp.get(),nquad1);
533 
534  // Inner product with respect to 'b' direction
535  for (mode=i=0; i < nmodes0; ++i)
536  {
537  Blas::Dgemv('T',nquad1,nmodes1-i,1.0, base1.get()+mode*nquad1,
538  nquad1,wsp.get()+i*nquad1,1, 0.0,
539  outarray.get() + mode,1);
540  mode += nmodes1 - i;
541  }
542 
543  // fix for modified basis by splitting top vertex mode
545  {
546  outarray[1] += Blas::Ddot(nquad1,base1.get()+nquad1,1,
547  wsp.get()+nquad1,1);
548  }
549  }
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:434
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:191
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 551 of file StdTriExp.cpp.

References v_IProductWRTDerivBase_SumFac().

555  {
556  StdTriExp::v_IProductWRTDerivBase_SumFac(dir,inarray,outarray);
557  }
virtual void v_IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdTriExp.cpp:593
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 559 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.

563  {
564  int nq = GetTotPoints();
566 
567  switch(dir)
568  {
569  case 0:
570  {
571  mtype = eIProductWRTDerivBase0;
572  break;
573  }
574  case 1:
575  {
576  mtype = eIProductWRTDerivBase1;
577  break;
578  }
579  default:
580  {
581  ASSERTL1(false,"input dir is out of range");
582  break;
583  }
584  }
585 
586  StdMatrixKey iprodmatkey(mtype,DetShapeType(),*this);
587  DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
588 
589  Blas::Dgemv('N',m_ncoeffs,nq,1.0,iprodmat->GetPtr().get(),
590  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
591  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:470
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:700
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:141
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
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 593 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().

597  {
598  int i;
599  int nquad0 = m_base[0]->GetNumPoints();
600  int nquad1 = m_base[1]->GetNumPoints();
601  int nqtot = nquad0*nquad1;
602  int nmodes0 = m_base[0]->GetNumModes();
603  int wspsize = max(max(nqtot,m_ncoeffs),nquad1*nmodes0);
604 
605  Array<OneD, NekDouble> gfac0(2*wspsize);
606  Array<OneD, NekDouble> tmp0 (gfac0+wspsize);
607 
608  const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
609 
610  // set up geometric factor: 2/(1-z1)
611  for(i = 0; i < nquad1; ++i)
612  {
613  gfac0[i] = 2.0/(1-z1[i]);
614  }
615 
616  for(i = 0; i < nquad1; ++i)
617  {
618  Vmath::Smul(nquad0,gfac0[i],&inarray[0]+i*nquad0,1,
619  &tmp0[0]+i*nquad0,1);
620  }
621 
622  MultiplyByQuadratureMetric(tmp0,tmp0);
623 
624  switch(dir)
625  {
626  case 0:
627  {
628  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
629  m_base[1]->GetBdata(),
630  tmp0,outarray,gfac0);
631  break;
632  }
633  case 1:
634  {
635  Array<OneD, NekDouble> tmp3(m_ncoeffs);
636  const Array<OneD, const NekDouble>& z0 = m_base[0]->GetZ();
637 
638  for (i = 0; i < nquad0; ++i)
639  {
640  gfac0[i] = 0.5*(1+z0[i]);
641  }
642 
643  for (i = 0; i < nquad1; ++i)
644  {
645  Vmath::Vmul(nquad0,&gfac0[0],1,&tmp0[0]+i*nquad0,1,
646  &tmp0[0]+i*nquad0,1);
647  }
648 
649  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
650  m_base[1]->GetBdata(),
651  tmp0,tmp3,gfac0);
652 
653  MultiplyByQuadratureMetric(inarray,tmp0);
654  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
655  m_base[1]->GetDbdata(),
656  tmp0,outarray,gfac0);
657  Vmath::Vadd(m_ncoeffs,&tmp3[0],1,&outarray[0],1,
658  &outarray[0],1);
659  break;
660  }
661  default:
662  {
663  ASSERTL1(false, "input dir is out of range");
664  break;
665  }
666  }
667  }
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:942
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
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:191
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:285
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:169
bool Nektar::StdRegions::StdTriExp::v_IsBoundaryInteriorExpansion ( )
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 849 of file StdTriExp.cpp.

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

850  {
851  return m_base[0]->GetBasisType() == LibUtilities::eModified_A &&
852  m_base[1]->GetBasisType() == LibUtilities::eModified_B;
853  }
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 1351 of file StdTriExp.cpp.

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

1355  {
1356  StdTriExp::v_LaplacianMatrixOp_MatFree(inarray,outarray,mkey);
1357  }
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 1359 of file StdTriExp.cpp.

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

1365  {
1367  k1,k2,inarray,outarray,mkey);
1368  }
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 673 of file StdTriExp.cpp.

References Nektar::NekConstants::kNekZeroTol.

675  {
676 
677  // set up local coordinate system
678  if (fabs(xi[1]-1.0) < NekConstants::kNekZeroTol)
679  {
680  eta[0] = -1.0;
681  eta[1] = 1.0;
682  }
683  else
684  {
685  eta[0] = 2*(1+xi[0])/(1-xi[1])-1.0;
686  eta[1] = xi[1];
687  }
688  }
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 1343 of file StdTriExp.cpp.

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

1347  {
1348  StdExpansion::MassMatrixOp_MatFree(inarray,outarray,mkey);
1349  }
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 1563 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().

1566  {
1567  int i;
1568  int nquad0 = m_base[0]->GetNumPoints();
1569  int nquad1 = m_base[1]->GetNumPoints();
1570 
1571  const Array<OneD, const NekDouble>& w0 = m_base[0]->GetW();
1572  const Array<OneD, const NekDouble>& w1 = m_base[1]->GetW();
1573  const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
1574 
1575  // multiply by integration constants
1576  for(i = 0; i < nquad1; ++i)
1577  {
1578  Vmath::Vmul(nquad0,inarray.get()+i*nquad0,1,
1579  w0.get(),1, outarray.get()+i*nquad0,1);
1580  }
1581 
1582  switch(m_base[1]->GetPointsType())
1583  {
1584  // Legendre inner product
1587  for(i = 0; i < nquad1; ++i)
1588  {
1589  Blas::Dscal(nquad0,0.5*(1-z1[i])*w1[i],
1590  outarray.get()+i*nquad0,1);
1591  }
1592  break;
1593 
1594  // (1,0) Jacobi Inner product
1596  for(i = 0; i < nquad1; ++i)
1597  {
1598  Blas::Dscal(nquad0,0.5*w1[i],outarray.get()+i*nquad0,1);
1599  }
1600  break;
1601 
1602  default:
1603  ASSERTL0(false, "Unsupported quadrature points type.");
1604  break;
1605  }
1606  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:57
1D Evenly-spaced points using Lagrange polynomial
Definition: PointsType.h:63
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:50
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:169
int Nektar::StdRegions::StdTriExp::v_NumBndryCoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp.

Definition at line 754 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().

755  {
757  "BasisType is not a boundary interior form");
759  "BasisType is not a boundary interior form");
760 
761  return 3 + (GetBasisNumModes(0)-2) + 2*(GetBasisNumModes(1)-2);
762  }
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:191
int Nektar::StdRegions::StdTriExp::v_NumDGBndryCoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 764 of file StdTriExp.cpp.

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

765  {
767  "BasisType is not a boundary interior form");
769  "BasisType is not a boundary interior form");
770 
771  return GetBasisNumModes(0) + 2*GetBasisNumModes(1);
772  }
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:191
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 128 of file StdTriExp.cpp.

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

Referenced by v_PhysDeriv(), and v_StdPhysDeriv().

133  {
134  int i;
135  int nquad0 = m_base[0]->GetNumPoints();
136  int nquad1 = m_base[1]->GetNumPoints();
137  Array<OneD, NekDouble> wsp(nquad0*nquad1);
138 
139  const Array<OneD, const NekDouble>& z0 = m_base[0]->GetZ();
140  const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
141 
142  // set up geometric factor: 2/(1-z1)
143  for (i = 0; i < nquad1; ++i)
144  {
145  wsp[i] = 2.0/(1-z1[i]);
146  }
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  for (i = 0; i < nquad0; ++i)
162  {
163  wsp[i] = 0.5*(1+z0[i]);
164  }
165 
166  for (i = 0; i < nquad1; ++i)
167  {
168  Vmath::Vvtvp(nquad0,&wsp[0],1,&out_d0[0]+i*nquad0,
169  1,&out_d1[0]+i*nquad0,
170  1,&out_d1[0]+i*nquad0,1);
171  }
172  }
173  }
174  else if (out_d1.num_elements() > 0)
175  {
176  Array<OneD, NekDouble> diff0(nquad0*nquad1);
177  PhysTensorDeriv(inarray, diff0, out_d1);
178 
179  for (i = 0; i < nquad1; ++i)
180  {
181  Blas::Dscal(nquad0,wsp[i],&diff0[0]+i*nquad0,1);
182  }
183 
184  for (i = 0; i < nquad0; ++i)
185  {
186  wsp[i] = 0.5*(1+z0[i]);
187  }
188 
189  for (i = 0; i < nquad1; ++i)
190  {
191  Vmath::Vvtvp(nquad0,&wsp[0],1,&diff0[0]+i*nquad0,
192  1,&out_d1[0]+i*nquad0,
193  1,&out_d1[0]+i*nquad0,1);
194  }
195  }
196  }
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:428
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...
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 198 of file StdTriExp.cpp.

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

202  {
203  switch(dir)
204  {
205  case 0:
206  {
207  v_PhysDeriv(inarray, outarray, NullNekDouble1DArray);
208  break;
209  }
210  case 1:
211  {
212  v_PhysDeriv(inarray, NullNekDouble1DArray, outarray);
213  break;
214  }
215  default:
216  {
217  ASSERTL1(false,"input dir is out of range");
218  break;
219  }
220  }
221  }
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:128
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
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 1480 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.

1484  {
1485  int n_coeffs = inarray.num_elements();
1486  int nquad0 = m_base[0]->GetNumPoints();
1487  int nquad1 = m_base[1]->GetNumPoints();
1488  Array<OneD, NekDouble> coeff(n_coeffs);
1489  Array<OneD, NekDouble> coeff_tmp(n_coeffs,0.0);
1490  Array<OneD, NekDouble> tmp;
1491  Array<OneD, NekDouble> tmp2;
1492  int nqtot = nquad0*nquad1;
1493  Array<OneD, NekDouble> phys_tmp(nqtot,0.0);
1494 
1495  int nmodes0 = m_base[0]->GetNumModes();
1496  int nmodes1 = m_base[1]->GetNumModes();
1497  int numMin2 = nmodes0;
1498  int i;
1499 
1500  const LibUtilities::PointsKey Pkey0(
1502  const LibUtilities::PointsKey Pkey1(
1504 
1505  LibUtilities::BasisKey b0(m_base[0]->GetBasisType(),nmodes0,Pkey0);
1506  LibUtilities::BasisKey b1(m_base[1]->GetBasisType(),nmodes1,Pkey1);
1507 
1508  LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A,nmodes0,Pkey0);
1509  LibUtilities::BasisKey bortho1(LibUtilities::eOrtho_B,nmodes1,Pkey1);
1510 
1511  StdRegions::StdTriExpSharedPtr m_OrthoTriExp;
1513 
1515  ::AllocateSharedPtr(b0, b1);
1516  m_OrthoTriExp = MemoryManager<StdRegions::StdTriExp>
1517  ::AllocateSharedPtr(bortho0, bortho1);
1518 
1519  m_TriExp ->BwdTrans(inarray,phys_tmp);
1520  m_OrthoTriExp->FwdTrans(phys_tmp, coeff);
1521 
1522  for (i = 0; i < n_coeffs; i++)
1523  {
1524  if (i == numMin)
1525  {
1526  coeff[i] = 0.0;
1527  numMin += numMin2 - 1;
1528  numMin2 -= 1.0;
1529  }
1530  }
1531 
1532  m_OrthoTriExp->BwdTrans(coeff,phys_tmp);
1533  m_TriExp ->FwdTrans(phys_tmp, outarray);
1534 
1535  }
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:50
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 223 of file StdTriExp.cpp.

References v_PhysDeriv().

228  {
229  StdTriExp::v_PhysDeriv(inarray, out_d0, out_d1);
230  }
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:128
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 232 of file StdTriExp.cpp.

References v_PhysDeriv().

236  {
237  StdTriExp::v_PhysDeriv(dir,inarray,outarray);
238  }
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:128
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 1388 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().

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

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

1375  {
1376  StdExpansion::WeakDerivMatrixOp_MatFree(i,inarray,outarray,mkey);
1377  }
void WeakDerivMatrixOp_MatFree(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)