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

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

Definition at line 52 of file StdTriExp.cpp.

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

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

Definition at line 69 of file StdTriExp.cpp.

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

Definition at line 75 of file StdTriExp.cpp.

76  {
77  }

Member Function Documentation

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

Backward tranform for triangular elements.

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

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::StdRegions::StdNodalTriExp.

Definition at line 252 of file StdTriExp.cpp.

References v_BwdTrans_SumFac().

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

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

280  {
281  int i;
282  int mode;
283  int nquad0 = m_base[0]->GetNumPoints();
284  int nquad1 = m_base[1]->GetNumPoints();
285  int nmodes0 = m_base[0]->GetNumModes();
286  int nmodes1 = m_base[1]->GetNumModes();
287 
288  ASSERTL1(wsp.num_elements() >= nquad0*nmodes1,
289  "Workspace size is not sufficient");
292  "Basis[1] is not of general tensor type");
293 
294  for (i = mode = 0; i < nmodes0; ++i)
295  {
296  Blas::Dgemv('N', nquad1,nmodes1-i,1.0,base1.get()+mode*nquad1,
297  nquad1,&inarray[0]+mode,1,0.0,&wsp[0]+i*nquad1,1);
298  mode += nmodes1-i;
299  }
300 
301  // fix for modified basis by splitting top vertex mode
303  {
304  Blas::Daxpy(nquad1,inarray[1],base1.get()+nquad1,1,
305  &wsp[0]+nquad1,1);
306  }
307 
308  Blas::Dgemm('N','T', nquad0,nquad1,nmodes0,1.0, base0.get(),nquad0,
309  &wsp[0], nquad1,0.0, &outarray[0], nquad0);
310  }
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 804 of file StdTriExp.cpp.

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

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

References v_GenMatrix().

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 857 of file StdTriExp.cpp.

References ASSERTL2.

858  {
859  ASSERTL2(edge >= 0 && edge <= 2, "edge id is out of range");
860 
861  return edge == 0 ? 0 : 1;
862  }
#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 864 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().

866  {
867  ASSERTL2(i >= 0 && i <= 2, "edge id is out of range");
868 
869  if (i == 0)
870  {
871  return GetBasis(0)->GetBasisKey();
872  }
873  else
874  {
875  switch(m_base[1]->GetBasisType())
876  {
878  {
879  switch(m_base[1]->GetPointsType())
880  {
882  {
883  LibUtilities::PointsKey pkey(
884  m_base[1]->GetBasisKey().GetPointsKey().
885  GetNumPoints()+1,
887  return LibUtilities::BasisKey(
889  m_base[1]->GetNumModes(),pkey);
890  break;
891  }
892 
893  default:
894  ASSERTL0(false,"unexpected points distribution");
895  break;
896  }
897  }
898  default:
899  ASSERTL0(false,"Information not available to set edge key");
900  break;
901  }
902  }
904  }
#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 751 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 692 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().

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

315  {
316  v_IProductWRTBase(inarray,outarray);
317 
318  // get Mass matrix inverse
319  StdMatrixKey masskey(eInvMass,DetShapeType(),*this);
320  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
321 
322  // copy inarray in case inarray == outarray
323  NekVector<NekDouble> in(m_ncoeffs,outarray,eCopy);
324  NekVector<NekDouble> out(m_ncoeffs,outarray,eWrapper);
325 
326  out = (*matsys)*in;
327  }
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:466
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 330 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().

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

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

1543  {
1545 
1546  if(inarray.get() == outarray.get())
1547  {
1548  Array<OneD,NekDouble> tmp(m_ncoeffs);
1549  Vmath::Vcopy(m_ncoeffs,inarray.get(),1,tmp.get(),1);
1550 
1551  Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, 1.0, mat->GetPtr().get(),
1552  m_ncoeffs, tmp.get(), 1, 0.0, outarray.get(), 1);
1553  }
1554  else
1555  {
1556  Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, 1.0, mat->GetPtr().get(),
1557  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
1558  }
1559  }
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 1258 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().

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

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

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

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

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

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

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

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

777  {
778  ASSERTL2(i >= 0 && i <= 2, "edge id is out of range");
779 
780  if (i == 0)
781  {
782  return GetBasisNumModes(0);
783  }
784  else
785  {
786  return GetBasisNumModes(1);
787  }
788  }
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 790 of file StdTriExp.cpp.

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

791  {
792  ASSERTL2((i >= 0)&&(i <= 2),"edge id is out of range");
793 
794  if (i == 0)
795  {
796  return GetNumPoints(0);
797  }
798  else
799  {
800  return GetNumPoints(1);
801  }
802  }
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 912 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().

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

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

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 741 of file StdTriExp.cpp.

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1610 of file StdTriExp.cpp.

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

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

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

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

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

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

Integrates the specified function over the domain.

See also
StdRegions::StdExpansion::Integral.

Reimplemented from Nektar::StdRegions::StdExpansion.

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

Definition at line 82 of file StdTriExp.cpp.

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

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

References v_IProductWRTBase_SumFac().

Referenced by v_FwdTrans(), and v_FwdTrans_BndConstrained().

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

476  {
477  int nq = GetTotPoints();
478  StdMatrixKey iprodmatkey(eIProductWRTBase,DetShapeType(),*this);
479  DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
480 
481  Blas::Dgemv('N',m_ncoeffs,nq,1.0,iprodmat->GetPtr().get(),
482  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
483  }
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 485 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().

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

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

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

References v_IProductWRTDerivBase_SumFac().

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

565  {
566  int nq = GetTotPoints();
568 
569  switch(dir)
570  {
571  case 0:
572  {
573  mtype = eIProductWRTDerivBase0;
574  break;
575  }
576  case 1:
577  {
578  mtype = eIProductWRTDerivBase1;
579  break;
580  }
581  default:
582  {
583  ASSERTL1(false,"input dir is out of range");
584  break;
585  }
586  }
587 
588  StdMatrixKey iprodmatkey(mtype,DetShapeType(),*this);
589  DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
590 
591  Blas::Dgemv('N',m_ncoeffs,nq,1.0,iprodmat->GetPtr().get(),
592  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
593  }
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 595 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().

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

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

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

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

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

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

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

References Nektar::NekConstants::kNekZeroTol.

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

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

1349  {
1350  StdExpansion::MassMatrixOp_MatFree(inarray,outarray,mkey);
1351  }
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 1565 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().

1568  {
1569  int i;
1570  int nquad0 = m_base[0]->GetNumPoints();
1571  int nquad1 = m_base[1]->GetNumPoints();
1572 
1573  const Array<OneD, const NekDouble>& w0 = m_base[0]->GetW();
1574  const Array<OneD, const NekDouble>& w1 = m_base[1]->GetW();
1575  const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
1576 
1577  // multiply by integration constants
1578  for(i = 0; i < nquad1; ++i)
1579  {
1580  Vmath::Vmul(nquad0,inarray.get()+i*nquad0,1,
1581  w0.get(),1, outarray.get()+i*nquad0,1);
1582  }
1583 
1584  switch(m_base[1]->GetPointsType())
1585  {
1586  // Legendre inner product
1589  for(i = 0; i < nquad1; ++i)
1590  {
1591  Blas::Dscal(nquad0,0.5*(1-z1[i])*w1[i],
1592  outarray.get()+i*nquad0,1);
1593  }
1594  break;
1595 
1596  // (1,0) Jacobi Inner product
1598  for(i = 0; i < nquad1; ++i)
1599  {
1600  Blas::Dscal(nquad0,0.5*w1[i],outarray.get()+i*nquad0,1);
1601  }
1602  break;
1603 
1604  default:
1605  ASSERTL0(false, "Unsupported quadrature points type.");
1606  break;
1607  }
1608  }
#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 756 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().

757  {
759  "BasisType is not a boundary interior form");
761  "BasisType is not a boundary interior form");
762 
763  return 3 + (GetBasisNumModes(0)-2) + 2*(GetBasisNumModes(1)-2);
764  }
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 766 of file StdTriExp.cpp.

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

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

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

Referenced by v_PhysDeriv(), and v_StdPhysDeriv().

135  {
136  int i;
137  int nquad0 = m_base[0]->GetNumPoints();
138  int nquad1 = m_base[1]->GetNumPoints();
139  Array<OneD, NekDouble> wsp(nquad0*nquad1);
140 
141  const Array<OneD, const NekDouble>& z0 = m_base[0]->GetZ();
142  const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
143 
144  // set up geometric factor: 2/(1-z1)
145  for (i = 0; i < nquad1; ++i)
146  {
147  wsp[i] = 2.0/(1-z1[i]);
148  }
149 
150  if (out_d0.num_elements() > 0)
151  {
152  PhysTensorDeriv(inarray, out_d0, out_d1);
153 
154  for (i = 0; i < nquad1; ++i)
155  {
156  Blas::Dscal(nquad0,wsp[i],&out_d0[0]+i*nquad0,1);
157  }
158 
159  // if no d1 required do not need to calculate both deriv
160  if (out_d1.num_elements() > 0)
161  {
162  // set up geometric factor: (1_z0)/(1-z1)
163  for (i = 0; i < nquad0; ++i)
164  {
165  wsp[i] = 0.5*(1+z0[i]);
166  }
167 
168  for (i = 0; i < nquad1; ++i)
169  {
170  Vmath::Vvtvp(nquad0,&wsp[0],1,&out_d0[0]+i*nquad0,
171  1,&out_d1[0]+i*nquad0,
172  1,&out_d1[0]+i*nquad0,1);
173  }
174  }
175  }
176  else if (out_d1.num_elements() > 0)
177  {
178  Array<OneD, NekDouble> diff0(nquad0*nquad1);
179  PhysTensorDeriv(inarray, diff0, out_d1);
180 
181  for (i = 0; i < nquad1; ++i)
182  {
183  Blas::Dscal(nquad0,wsp[i],&diff0[0]+i*nquad0,1);
184  }
185 
186  for (i = 0; i < nquad0; ++i)
187  {
188  wsp[i] = 0.5*(1+z0[i]);
189  }
190 
191  for (i = 0; i < nquad1; ++i)
192  {
193  Vmath::Vvtvp(nquad0,&wsp[0],1,&diff0[0]+i*nquad0,
194  1,&out_d1[0]+i*nquad0,
195  1,&out_d1[0]+i*nquad0,1);
196  }
197  }
198  }
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 200 of file StdTriExp.cpp.

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

204  {
205  switch(dir)
206  {
207  case 0:
208  {
209  v_PhysDeriv(inarray, outarray, NullNekDouble1DArray);
210  break;
211  }
212  case 1:
213  {
214  v_PhysDeriv(inarray, NullNekDouble1DArray, outarray);
215  break;
216  }
217  default:
218  {
219  ASSERTL1(false,"input dir is out of range");
220  break;
221  }
222  }
223  }
static Array< OneD, NekDouble > NullNekDouble1DArray
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Calculate the derivative of the physical points.
Definition: StdTriExp.cpp:130
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp: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 1482 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.

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

References v_PhysDeriv().

230  {
231  StdTriExp::v_PhysDeriv(inarray, out_d0, out_d1);
232  }
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Calculate the derivative of the physical points.
Definition: StdTriExp.cpp:130
void Nektar::StdRegions::StdTriExp::v_StdPhysDeriv ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 234 of file StdTriExp.cpp.

References v_PhysDeriv().

238  {
239  StdTriExp::v_PhysDeriv(dir,inarray,outarray);
240  }
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Calculate the derivative of the physical points.
Definition: StdTriExp.cpp:130
void Nektar::StdRegions::StdTriExp::v_SVVLaplacianFilter ( Array< OneD, NekDouble > &  array,
const StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::TriExp.

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

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

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

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