Nektar++
Public Member Functions | Protected Member Functions | Private Types | List of all members
Nektar::StdRegions::StdQuadExp Class Reference

#include <StdQuadExp.h>

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

Public Member Functions

 StdQuadExp ()
 
 StdQuadExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb)
 Constructor using BasisKey class for quadrature points and order definition. More...
 
 StdQuadExp (const StdQuadExp &T)
 Copy Constructor. More...
 
 ~StdQuadExp ()
 Destructor. More...
 
- 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::BasisSharedPtrGetBasis (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
 
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, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
void DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
StdRegions::Orientation 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)
 
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 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, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual void v_DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual StdRegions::Orientation v_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)
 
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)
 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...
 
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

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)
 
virtual void v_FwdTrans_BndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &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_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Calculate the inner product of inarray with respect to the basis B=base0*base1 and put into outarray. More...
 
virtual void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
virtual void v_IProductWRTBase_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &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_SumFac (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_LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 
virtual void v_FillMode (const int mode, Array< OneD, NekDouble > &array)
 Fill outarray with mode mode of expansion. More...
 
virtual int v_GetNverts () const
 
virtual int v_GetNedges () const
 
virtual int v_GetEdgeNcoeffs (const int i) const
 
virtual int v_GetEdgeNumPoints (const int i) const
 
virtual int v_NumBndryCoeffs () const
 
virtual int v_NumDGBndryCoeffs () const
 
virtual LibUtilities::BasisType v_GetEdgeBasisType (const int i) const
 
virtual int v_DetCartesianDirOfEdge (const int edge)
 
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
virtual const LibUtilities::BasisKey v_DetEdgeBasisKey (const int i) const
 
virtual LibUtilities::ShapeType v_DetShapeType () const
 
virtual bool v_IsBoundaryInteriorExpansion ()
 
virtual void v_GetCoords (Array< OneD, NekDouble > &coords_0, Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2)
 
virtual void v_GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
virtual void v_GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
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)
 
void v_GetEdgeToElementMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
 
virtual DNekMatSharedPtr v_GenMatrix (const StdMatrixKey &mkey)
 
virtual DNekMatSharedPtr v_CreateStdMatrix (const StdMatrixKey &mkey)
 
virtual void v_MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
virtual void v_ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void v_GeneralMatrixOp_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
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)
 
- 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)
 

Private Types

typedef boost::shared_ptr< StdExpansion1DStdExpansion1DSharedPtr
 

Additional Inherited Members

- Protected Attributes inherited from Nektar::StdRegions::StdExpansion
Array< OneD, LibUtilities::BasisSharedPtrm_base
 
int m_elmt_id
 
int m_ncoeffs
 
LibUtilities::NekManager< StdMatrixKey, DNekMat, StdMatrixKey::opLessm_stdMatrixManager
 
LibUtilities::NekManager< StdMatrixKey, DNekBlkMat, StdMatrixKey::opLessm_stdStaticCondMatrixManager
 
LibUtilities::NekManager< IndexMapKey, IndexMapValues, IndexMapKey::opLessm_IndexMapManager
 

Detailed Description

Definition at line 51 of file StdQuadExp.h.

Member Typedef Documentation

Definition at line 54 of file StdQuadExp.h.

Constructor & Destructor Documentation

Nektar::StdRegions::StdQuadExp::StdQuadExp ( )

Definition at line 46 of file StdQuadExp.cpp.

47  {
48  }
Nektar::StdRegions::StdQuadExp::StdQuadExp ( const LibUtilities::BasisKey Ba,
const LibUtilities::BasisKey Bb 
)

Constructor using BasisKey class for quadrature points and order definition.

Definition at line 53 of file StdQuadExp.cpp.

54  :
55  StdExpansion (Ba.GetNumModes()*Bb.GetNumModes(),2,Ba,Bb),
56  StdExpansion2D(Ba.GetNumModes()*Bb.GetNumModes(),Ba,Bb)
57  {
58  }
StdExpansion()
Default Constructor.
Nektar::StdRegions::StdQuadExp::StdQuadExp ( const StdQuadExp T)

Copy Constructor.

Definition at line 61 of file StdQuadExp.cpp.

61  :
62  StdExpansion(T),
64  {
65  }
StdExpansion()
Default Constructor.
Nektar::StdRegions::StdQuadExp::~StdQuadExp ( )

Destructor.

Definition at line 68 of file StdQuadExp.cpp.

69  {
70  }

Member Function Documentation

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

Implements Nektar::StdRegions::StdExpansion.

Definition at line 150 of file StdQuadExp.cpp.

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

152  {
153  if(m_base[0]->Collocation() && m_base[1]->Collocation())
154  {
156  inarray, 1, outarray, 1);
157  }
158  else
159  {
160  StdQuadExp::v_BwdTrans_SumFac(inarray,outarray);
161  }
162  }
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:229
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdQuadExp.cpp:164
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1038
void Nektar::StdRegions::StdQuadExp::v_BwdTrans_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 164 of file StdQuadExp.cpp.

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

Referenced by v_BwdTrans().

167  {
168  Array<OneD, NekDouble> wsp(m_base[0]->GetNumPoints()*
169  m_base[1]->GetNumModes());
170 
171  BwdTrans_SumFacKernel(m_base[0]->GetBdata(),
172  m_base[1]->GetBdata(),
173  inarray,outarray,wsp,true,true);
174  }
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::StdQuadExp::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 186 of file StdQuadExp.cpp.

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

194  {
195  int nquad0 = m_base[0]->GetNumPoints();
196  int nquad1 = m_base[1]->GetNumPoints();
197  int nmodes0 = m_base[0]->GetNumModes();
198  int nmodes1 = m_base[1]->GetNumModes();
199 
200  bool colldir0 = doCheckCollDir0?(m_base[0]->Collocation()):false;
201  bool colldir1 = doCheckCollDir1?(m_base[1]->Collocation()):false;
202 
203  if(colldir0 && colldir1)
204  {
205  Vmath::Vcopy(m_ncoeffs,inarray.get(),1,outarray.get(),1);
206  }
207  else if(colldir0)
208  {
209  Blas::Dgemm('N','T', nquad0, nquad1,nmodes1, 1.0, &inarray[0], nquad0,
210  base1.get(), nquad1, 0.0, &outarray[0], nquad0);
211  }
212  else if(colldir1)
213  {
214  Blas::Dgemm('N','N', nquad0,nmodes1,nmodes0,1.0, base0.get(),
215  nquad0, &inarray[0], nmodes0,0.0,&outarray[0], nquad0);
216  }
217  else
218  {
219  ASSERTL1(wsp.num_elements()>=nquad0*nmodes1,"Workspace size is not sufficient");
220 
221  // Those two calls correpsond to the operation
222  // out = B0*in*Transpose(B1);
223  Blas::Dgemm('N','N', nquad0,nmodes1,nmodes0,1.0, base0.get(),
224  nquad0, &inarray[0], nmodes0,0.0,&wsp[0], nquad0);
225  Blas::Dgemm('N','T', nquad0, nquad1,nmodes1, 1.0, &wsp[0], nquad0,
226  base1.get(), nquad1, 0.0, &outarray[0], nquad0);
227  }
228  }
#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 Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1038
int Nektar::StdRegions::StdQuadExp::v_CalcNumberOfCoefficients ( const std::vector< unsigned int > &  nummodes,
int &  modes_offset 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 731 of file StdQuadExp.cpp.

734  {
735  int nmodes = nummodes[modes_offset]*nummodes[modes_offset+1];
736  modes_offset += 2;
737 
738  return nmodes;
739  }
DNekMatSharedPtr Nektar::StdRegions::StdQuadExp::v_CreateStdMatrix ( const StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 1444 of file StdQuadExp.cpp.

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

1445  {
1446  return GenMatrix(mkey);
1447  }
DNekMatSharedPtr GenMatrix(const StdMatrixKey &mkey)
int Nektar::StdRegions::StdQuadExp::v_DetCartesianDirOfEdge ( const int  edge)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 668 of file StdQuadExp.cpp.

References ASSERTL2.

669  {
670  ASSERTL2((edge >= 0)&&(edge <= 3),"edge id is out of range");
671 
672  if((edge == 0)||(edge == 2))
673  {
674  return 0;
675  }
676  else
677  {
678  return 1;
679  }
680  }
#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::StdQuadExp::v_DetEdgeBasisKey ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 682 of file StdQuadExp.cpp.

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

683  {
684  ASSERTL2((i >= 0)&&(i <= 3),"edge id is out of range");
685 
686  if((i == 0)||(i == 2))
687  {
688  return GetBasis(0)->GetBasisKey();
689  }
690  else
691  {
692  return GetBasis(1)->GetBasisKey();
693  }
694 
695  }
const LibUtilities::BasisSharedPtr & GetBasis(int dir) const
This function gets the shared point to basis in the dir direction.
Definition: StdExpansion.h:118
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:213
LibUtilities::ShapeType Nektar::StdRegions::StdQuadExp::v_DetShapeType ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 697 of file StdQuadExp.cpp.

References Nektar::LibUtilities::eQuadrilateral.

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

Fill outarray with mode mode of expansion.

Note for quadrilateral expansions _base0 modes run fastest

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 580 of file StdQuadExp.cpp.

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

582  {
583  int i;
584  int nquad0 = m_base[0]->GetNumPoints();
585  int nquad1 = m_base[1]->GetNumPoints();
586  Array<OneD, const NekDouble> base0 = m_base[0]->GetBdata();
587  Array<OneD, const NekDouble> base1 = m_base[1]->GetBdata();
588  int btmp0 = m_base[0]->GetNumModes();
589  int mode0 = mode%btmp0;
590  int mode1 = mode/btmp0;
591 
592 
593  ASSERTL2(mode1 == (int)floor((1.0*mode)/btmp0),
594  "Integer Truncation not Equiv to Floor");
595 
596  ASSERTL2(m_ncoeffs <= mode,
597  "calling argument mode is larger than total expansion order");
598 
599  for(i = 0; i < nquad1; ++i)
600  {
601  Vmath::Vcopy(nquad0,(NekDouble *)(base0.get() + mode0*nquad0),
602  1, &outarray[0]+i*nquad0,1);
603  }
604 
605  for(i = 0; i < nquad0; ++i)
606  {
607  Vmath::Vmul(nquad1,(NekDouble *)(base1.get() + mode1*nquad1),1,
608  &outarray[0]+i,nquad0,&outarray[0]+i,nquad0);
609  }
610  }
double NekDouble
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp: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:1038
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::StdQuadExp::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::QuadExp.

Definition at line 231 of file StdQuadExp.cpp.

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

233  {
234  if((m_base[0]->Collocation())&&(m_base[1]->Collocation()))
235  {
236  Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
237  }
238  else
239  {
240  StdQuadExp::v_IProductWRTBase(inarray,outarray);
241 
242  // get Mass matrix inverse
243  StdMatrixKey masskey(eInvMass,DetShapeType(),*this);
244  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
245 
246  // copy inarray in case inarray == outarray
247  NekVector<NekDouble> in(m_ncoeffs,outarray,eCopy);
248  NekVector<NekDouble> out(m_ncoeffs,outarray,eWrapper);
249 
250  out = (*matsys)*in;
251  }
252  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:454
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:684
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*base1 and put into outarray...
Definition: StdQuadExp.cpp:379
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1038
void Nektar::StdRegions::StdQuadExp::v_FwdTrans_BndConstrained ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 254 of file StdQuadExp.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::StdRegions::eForwards, Nektar::StdRegions::eMass, Nektar::StdRegions::StdExpansion::GetEdgeToElementMap(), Nektar::StdRegions::StdExpansion::GetInteriorMap(), Nektar::StdRegions::StdExpansion::IProductWRTBase(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::StdRegions::StdExpansion::m_stdStaticCondMatrixManager, Nektar::StdRegions::StdExpansion::MassMatrixOp(), Nektar::StdRegions::StdExpansion::NumBndryCoeffs(), sign, Vmath::Vcopy(), and Vmath::Vsub().

257  {
258  if((m_base[0]->Collocation())&&(m_base[1]->Collocation()))
259  {
260  Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
261  }
262  else
263  {
264  int i,j;
265  int npoints[2] = {m_base[0]->GetNumPoints(),
266  m_base[1]->GetNumPoints()};
267  int nmodes[2] = {m_base[0]->GetNumModes(),
268  m_base[1]->GetNumModes()};
269 
270  fill(outarray.get(), outarray.get()+m_ncoeffs, 0.0 );
271 
272  Array<OneD, NekDouble> physEdge[4];
273  Array<OneD, NekDouble> coeffEdge[4];
274  for(i = 0; i < 4; i++)
275  {
276  physEdge[i] = Array<OneD, NekDouble>(npoints[i%2]);
277  coeffEdge[i] = Array<OneD, NekDouble>(nmodes[i%2]);
278  }
279 
280  for(i = 0; i < npoints[0]; i++)
281  {
282  physEdge[0][i] = inarray[i];
283  physEdge[2][i] = inarray[npoints[0]*npoints[1]-1-i];
284  }
285 
286  for(i = 0; i < npoints[1]; i++)
287  {
288  physEdge[1][i] = inarray[npoints[0]-1+i*npoints[0]];
289  physEdge[3][i] = inarray[(npoints[1]-1)*npoints[0]-i*npoints[0]];
290  }
291 
294 
295  Array<OneD, unsigned int> mapArray;
296  Array<OneD, int> signArray;
297  NekDouble sign;
298 
299  for(i = 0; i < 4; i++)
300  {
301  segexp[i%2]->FwdTrans_BndConstrained(physEdge[i],coeffEdge[i]);
302 
303  GetEdgeToElementMap(i,eForwards,mapArray,signArray);
304  for(j=0; j < nmodes[i%2]; j++)
305  {
306  sign = (NekDouble) signArray[j];
307  outarray[ mapArray[j] ] = sign * coeffEdge[i][j];
308  }
309  }
310 
311  Array<OneD, NekDouble> tmp0(m_ncoeffs);
312  Array<OneD, NekDouble> tmp1(m_ncoeffs);
313 
314  StdMatrixKey masskey(eMass,DetShapeType(),*this);
315  MassMatrixOp(outarray,tmp0,masskey);
316  IProductWRTBase(inarray,tmp1);
317 
318  Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
319 
320  // get Mass matrix inverse (only of interior DOF)
321  // use block (1,1) of the static condensed system
322  // note: this block alreay contains the inverse matrix
323  DNekMatSharedPtr matsys = (m_stdStaticCondMatrixManager[masskey])->GetBlock(1,1);
324 
325  int nBoundaryDofs = NumBndryCoeffs();
326  int nInteriorDofs = m_ncoeffs - nBoundaryDofs;
327 
328 
329  Array<OneD, NekDouble> rhs(nInteriorDofs);
330  Array<OneD, NekDouble> result(nInteriorDofs);
331 
332  GetInteriorMap(mapArray);
333 
334  for(i = 0; i < nInteriorDofs; i++)
335  {
336  rhs[i] = tmp1[ mapArray[i] ];
337  }
338 
339  Blas::Dgemv('N',nInteriorDofs,nInteriorDofs,1.0, &(matsys->GetPtr())[0],
340  nInteriorDofs,rhs.get(),1,0.0,result.get(),1);
341 
342  for(i = 0; i < nInteriorDofs; i++)
343  {
344  outarray[ mapArray[i] ] = result[i];
345  }
346  }
347 
348  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:454
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:931
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
void GetEdgeToElementMap(const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
Definition: StdExpansion.h:820
void IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
this function calculates the inner product of a given function f with the different modes of the expa...
Definition: StdExpansion.h:613
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
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:795
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 Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1038
void Nektar::StdRegions::StdQuadExp::v_GeneralMatrixOp_MatOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdMatrixKey mkey 
)
protected

Definition at line 1453 of file StdQuadExp.cpp.

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

1457  {
1458 
1460 
1461  if(inarray.get() == outarray.get())
1462  {
1463  Array<OneD,NekDouble> tmp(m_ncoeffs);
1464  Vmath::Vcopy(m_ncoeffs,inarray.get(),1,tmp.get(),1);
1465 
1466  Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, 1.0, mat->GetPtr().get(),
1467  m_ncoeffs, tmp.get(), 1, 0.0, outarray.get(), 1);
1468  }
1469  else
1470  {
1471  Blas::Dgemv('N', m_ncoeffs, m_ncoeffs, 1.0, mat->GetPtr().get(),
1472  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
1473  }
1474  }
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:1038
DNekMatSharedPtr Nektar::StdRegions::StdQuadExp::v_GenMatrix ( const StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 1319 of file StdQuadExp.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::LibUtilities::BasisManager(), Nektar::StdRegions::StdExpansion::CreateGeneralMatrix(), Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::StdRegions::eFactorGaussEdge, Nektar::LibUtilities::eFourier, Nektar::StdRegions::eFwdTrans, Nektar::LibUtilities::eGauss_Lagrange, Nektar::StdRegions::eGaussDG, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::StdRegions::eInvMass, Nektar::StdRegions::eIProductWRTBase, Nektar::StdRegions::eMass, Nektar::StdRegions::ePhysInterpToEquiSpaced, Nektar::StdRegions::StdExpansion::GetBasisNumModes(), Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdMatrixKey::GetConstFactors(), Nektar::StdRegions::StdMatrixKey::GetMatrixType(), Nektar::LibUtilities::StdQuadData::getNumberOfCoefficients(), Nektar::StdRegions::StdExpansion::GetStdMatrix(), Nektar::StdRegions::StdExpansion::LocCoordToLocCollapsed(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Vmath::Smul(), and Vmath::Vcopy().

1320  {
1321  int i;
1322  int order0 = GetBasisNumModes(0);
1323  int order1 = GetBasisNumModes(1);
1324  MatrixType mtype = mkey.GetMatrixType();
1325 
1326  DNekMatSharedPtr Mat;
1327 
1328  switch(mtype)
1329  {
1331  {
1332  int nq0 = m_base[0]->GetNumPoints();
1333  int nq1 = m_base[1]->GetNumPoints();
1334  int nq = max(nq0,nq1);
1335  int neq = LibUtilities::StdQuadData::
1336  getNumberOfCoefficients(nq, nq);
1337  Array<OneD, Array<OneD, NekDouble> > coords(neq);
1338  Array<OneD, NekDouble> coll (2);
1339  Array<OneD, DNekMatSharedPtr> I (2);
1340  Array<OneD, NekDouble> tmp (nq0);
1341 
1343  AllocateSharedPtr(neq, nq0 * nq1);
1344  int cnt = 0;
1345 
1346  for(int i = 0; i < nq; ++i)
1347  {
1348  for(int j = 0; j < nq; ++j,++cnt)
1349  {
1350  coords[cnt] = Array<OneD, NekDouble>(2);
1351  coords[cnt][0] = -1.0 + 2*j/(NekDouble)(nq-1);
1352  coords[cnt][1] = -1.0 + 2*i/(NekDouble)(nq-1);
1353  }
1354  }
1355 
1356  for(int i = 0; i < neq; ++i)
1357  {
1358  LocCoordToLocCollapsed(coords[i],coll);
1359 
1360  I[0] = m_base[0]->GetI(coll);
1361  I[1] = m_base[1]->GetI(coll+1);
1362 
1363  // interpolate first coordinate direction
1364  for (int j = 0; j < nq1; ++j)
1365  {
1366  NekDouble fac = (I[1]->GetPtr())[j];
1367  Vmath::Smul(nq0,fac,I[0]->GetPtr(),1,tmp,1);
1368 
1369  Vmath::Vcopy(nq0, &tmp[0], 1,
1370  Mat->GetRawPtr()+j*nq0*neq+i,neq);
1371  }
1372 
1373  }
1374  break;
1375  }
1376  case eMass:
1377  {
1379  // For Fourier basis set the imaginary component of mean mode
1380  // to have a unit diagonal component in mass matrix
1382  {
1383  for(i = 0; i < order1; ++i)
1384  {
1385  (*Mat)(order0*i+1,i*order0+1) = 1.0;
1386  }
1387  }
1388 
1390  {
1391  for(i = 0; i < order0; ++i)
1392  {
1393  (*Mat)(order0+i ,order0+i) = 1.0;
1394  }
1395  }
1396  break;
1397  }
1398  case eFwdTrans:
1399  {
1401  StdMatrixKey iprodkey(eIProductWRTBase,DetShapeType(),*this);
1402  DNekMat &Iprod = *GetStdMatrix(iprodkey);
1403  StdMatrixKey imasskey(eInvMass,DetShapeType(),*this);
1404  DNekMat &Imass = *GetStdMatrix(imasskey);
1405 
1406  (*Mat) = Imass*Iprod;
1407  break;
1408  }
1409  case eGaussDG:
1410  {
1411  ConstFactorMap factors = mkey.GetConstFactors();
1412 
1413  int edge = (int)factors[StdRegions::eFactorGaussEdge];
1414  int dir = (edge + 1) % 2;
1415  int nCoeffs = m_base[dir]->GetNumModes();
1416 
1417  const LibUtilities::PointsKey BS_p(
1419  const LibUtilities::BasisKey BS_k(
1420  LibUtilities::eGauss_Lagrange, nCoeffs, BS_p);
1421 
1422  Array<OneD, NekDouble> coords(1, 0.0);
1423  coords[0] = (edge == 0 || edge == 3) ? -1.0 : 1.0;
1424 
1427  DNekMatSharedPtr m_Ix = basis->GetI(coords);
1428 
1430  1.0, nCoeffs);
1431  Vmath::Vcopy(nCoeffs, m_Ix->GetPtr(), 1, Mat->GetPtr(), 1);
1432  break;
1433  }
1434  default:
1435  {
1437  break;
1438  }
1439  }
1440 
1441  return Mat;
1442  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:454
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
int GetBasisNumModes(const int dir) const
This function returns the number of expansion modes in the dir direction.
Definition: StdExpansion.h:178
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:54
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:248
Fourier Expansion .
Definition: BasisType.h:52
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:684
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:47
BasisManagerT & BasisManager(void)
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
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
NekMatrix< NekDouble, StandardMatrixTag > DNekMat
Definition: NekTypeDefs.hpp:52
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
int getNumberOfCoefficients(int Na, int Nb)
Definition: ShapeType.hpp:132
boost::shared_ptr< Basis > BasisSharedPtr
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1038
void Nektar::StdRegions::StdQuadExp::v_GetBoundaryMap ( Array< OneD, unsigned int > &  outarray)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 779 of file StdQuadExp.cpp.

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

780  {
781  int i;
782  int cnt=0;
783  int nummodes0, nummodes1;
784  int value1 = 0, value2 = 0;
785  if(outarray.num_elements()!=NumBndryCoeffs())
786  {
787  outarray = Array<OneD, unsigned int>(NumBndryCoeffs());
788  }
789 
790  nummodes0 = m_base[0]->GetNumModes();
791  nummodes1 = m_base[1]->GetNumModes();
792 
793  const LibUtilities::BasisType Btype0 = GetBasisType(0);
794  const LibUtilities::BasisType Btype1 = GetBasisType(1);
795 
796  switch(Btype1)
797  {
800  value1 = nummodes0;
801  break;
803  value1 = 2*nummodes0;
804  break;
805  default:
806  ASSERTL0(0,"Mapping array is not defined for this expansion");
807  break;
808  }
809 
810  for(i = 0; i < value1; i++)
811  {
812  outarray[i]=i;
813  }
814  cnt=value1;
815 
816  switch(Btype0)
817  {
820  value2 = value1+nummodes0-1;
821  break;
823  value2 = value1+1;
824  break;
825  default:
826  ASSERTL0(0,"Mapping array is not defined for this expansion");
827  break;
828  }
829 
830  for(i = 0; i < nummodes1-2; i++)
831  {
832  outarray[cnt++]=value1+i*nummodes0;
833  outarray[cnt++]=value2+i*nummodes0;
834  }
835 
836 
838  {
839  for(i = nummodes0*(nummodes1-1);i < GetNcoeffs(); i++)
840  {
841  outarray[cnt++] = i;
842  }
843  }
844  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
Principle Modified Functions .
Definition: BasisType.h:49
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:54
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
Lagrange for SEM basis .
Definition: BasisType.h:53
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdQuadExp::v_GetCoords ( Array< OneD, NekDouble > &  coords_0,
Array< OneD, NekDouble > &  coords_1,
Array< OneD, NekDouble > &  coords_2 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 758 of file StdQuadExp.cpp.

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

761  {
762  Array<OneD, const NekDouble> z0 = m_base[0]->GetZ();
763  Array<OneD, const NekDouble> z1 = m_base[1]->GetZ();
764  int nq0 = GetNumPoints(0);
765  int nq1 = GetNumPoints(1);
766  int i;
767 
768  for(i = 0; i < nq1; ++i)
769  {
770  Blas::Dcopy(nq0,z0.get(), 1,&coords_0[0] + i*nq0,1);
771  Vmath::Fill(nq0,z1[i],&coords_1[0] + i*nq0,1);
772  }
773  }
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::StdQuadExp::v_GetEdgeBasisType ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 654 of file StdQuadExp.cpp.

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

655  {
656  ASSERTL2((i >= 0)&&(i <= 3),"edge id is out of range");
657 
658  if((i == 0)||(i == 2))
659  {
660  return GetBasisType(0);
661  }
662  else
663  {
664  return GetBasisType(1);
665  }
666  }
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::StdQuadExp::v_GetEdgeInteriorMap ( const int  eid,
const Orientation  edgeOrient,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1007 of file StdQuadExp.cpp.

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

1011  {
1012  int i;
1013  const int nummodes0 = m_base[0]->GetNumModes();
1014  const int nummodes1 = m_base[1]->GetNumModes();
1015  const int nEdgeIntCoeffs = GetEdgeNcoeffs(eid)-2;
1016  const LibUtilities::BasisType bType = GetEdgeBasisType(eid);
1017 
1018  if(maparray.num_elements() != nEdgeIntCoeffs)
1019  {
1020  maparray = Array<OneD, unsigned int>(nEdgeIntCoeffs);
1021  }
1022 
1023  if(signarray.num_elements() != nEdgeIntCoeffs)
1024  {
1025  signarray = Array<OneD, int>(nEdgeIntCoeffs,1);
1026  }
1027  else
1028  {
1029  fill( signarray.get() , signarray.get()+nEdgeIntCoeffs, 1 );
1030  }
1031 
1032  if(bType == LibUtilities::eModified_A)
1033  {
1034  switch(eid)
1035  {
1036  case 0:
1037  {
1038  for(i = 0; i < nEdgeIntCoeffs; i++)
1039  {
1040  maparray[i] = i+2;
1041  }
1042 
1043  if(edgeOrient==eBackwards)
1044  {
1045  for(i = 1; i < nEdgeIntCoeffs; i+=2)
1046  {
1047  signarray[i] = -1;
1048  }
1049  }
1050  }
1051  break;
1052  case 1:
1053  {
1054  for(i = 0; i < nEdgeIntCoeffs; i++)
1055  {
1056  maparray[i] = (i+2)*nummodes0 + 1;
1057  }
1058 
1059  if(edgeOrient==eBackwards)
1060  {
1061  for(i = 1; i < nEdgeIntCoeffs; i+=2)
1062  {
1063  signarray[i] = -1;
1064  }
1065  }
1066  }
1067  break;
1068  case 2:
1069  {
1070  for(i = 0; i < nEdgeIntCoeffs; i++)
1071  {
1072  maparray[i] = nummodes0+i+2;
1073  }
1074 
1075  if(edgeOrient==eForwards)
1076  {
1077  for(i = 1; i < nEdgeIntCoeffs; i+=2)
1078  {
1079  signarray[i] = -1;
1080  }
1081  }
1082  }
1083  break;
1084  case 3:
1085  {
1086  for(i = 0; i < nEdgeIntCoeffs; i++)
1087  {
1088  maparray[i] = (i+2)*nummodes0;
1089  }
1090 
1091  if(edgeOrient==eForwards)
1092  {
1093  for(i = 1; i < nEdgeIntCoeffs; i+=2)
1094  {
1095  signarray[i] = -1;
1096  }
1097  }
1098  }
1099  break;
1100  default:
1101  ASSERTL0(false,"eid must be between 0 and 3");
1102  break;
1103  }
1104  }
1105  else if(bType == LibUtilities::eGLL_Lagrange)
1106  {
1107  switch(eid)
1108  {
1109  case 0:
1110  {
1111  for(i = 0; i < nEdgeIntCoeffs; i++)
1112  {
1113  maparray[i] = i+1;
1114  }
1115  }
1116  break;
1117  case 1:
1118  {
1119  for(i = 0; i < nEdgeIntCoeffs; i++)
1120  {
1121  maparray[i] = (i+2)*nummodes0 - 1;
1122  }
1123  }
1124  break;
1125  case 2:
1126  {
1127  for(i = 0; i < nEdgeIntCoeffs; i++)
1128  {
1129  maparray[i] = nummodes0*nummodes1 - 2 - i;
1130  }
1131  }
1132  break;
1133  case 3:
1134  {
1135  for(i = 0; i < nEdgeIntCoeffs; i++)
1136  {
1137  maparray[i] = nummodes0*(nummodes1-2-i);
1138  }
1139  }
1140  break;
1141  default:
1142  ASSERTL0(false,"eid must be between 0 and 3");
1143  break;
1144  }
1145  if(edgeOrient == eBackwards)
1146  {
1147  reverse( maparray.get() , maparray.get()+nEdgeIntCoeffs );
1148  }
1149  }
1150  else
1151  {
1152  ASSERTL0(false,"Mapping not defined for this type of basis");
1153  }
1154 
1155  }
#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:397
Lagrange for SEM basis .
Definition: BasisType.h:53
Array< OneD, LibUtilities::BasisSharedPtr > m_base
int Nektar::StdRegions::StdQuadExp::v_GetEdgeNcoeffs ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 626 of file StdQuadExp.cpp.

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

627  {
628  ASSERTL2((i >= 0)&&(i <= 3),"edge id is out of range");
629 
630  if((i == 0)||(i == 2))
631  {
632  return GetBasisNumModes(0);
633  }
634  else
635  {
636  return GetBasisNumModes(1);
637  }
638  }
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::StdQuadExp::v_GetEdgeNumPoints ( const int  i) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 640 of file StdQuadExp.cpp.

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

641  {
642  ASSERTL2((i >= 0)&&(i <= 3),"edge id is out of range");
643 
644  if((i == 0)||(i == 2))
645  {
646  return GetNumPoints(0);
647  }
648  else
649  {
650  return GetNumPoints(1);
651  }
652  }
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::StdQuadExp::v_GetEdgeToElementMap ( const int  eid,
const Orientation  edgeOrient,
Array< OneD, unsigned int > &  maparray,
Array< OneD, int > &  signarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1157 of file StdQuadExp.cpp.

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

1161  {
1162  int i;
1163  const int nummodes0 = m_base[0]->GetNumModes();
1164  const int nummodes1 = m_base[1]->GetNumModes();
1165  const int nEdgeCoeffs = GetEdgeNcoeffs(eid);
1166  const LibUtilities::BasisType bType = GetEdgeBasisType(eid);
1167 
1168  if(maparray.num_elements() != nEdgeCoeffs)
1169  {
1170  maparray = Array<OneD, unsigned int>(nEdgeCoeffs);
1171  }
1172 
1173  if(signarray.num_elements() != nEdgeCoeffs)
1174  {
1175  signarray = Array<OneD, int>(nEdgeCoeffs,1);
1176  }
1177  else
1178  {
1179  fill( signarray.get() , signarray.get()+nEdgeCoeffs, 1 );
1180  }
1181 
1182  if(bType == LibUtilities::eModified_A)
1183  {
1184  switch(eid)
1185  {
1186  case 0:
1187  {
1188  for(i = 0; i < nEdgeCoeffs; i++)
1189  {
1190  maparray[i] = i;
1191  }
1192 
1193  if(edgeOrient==eBackwards)
1194  {
1195  swap( maparray[0] , maparray[1] );
1196 
1197  for(i = 3; i < nEdgeCoeffs; i+=2)
1198  {
1199  signarray[i] = -1;
1200  }
1201  }
1202  }
1203  break;
1204  case 1:
1205  {
1206  for(i = 0; i < nEdgeCoeffs; i++)
1207  {
1208  maparray[i] = i*nummodes0 + 1;
1209  }
1210 
1211  if(edgeOrient==eBackwards)
1212  {
1213  swap( maparray[0] , maparray[1] );
1214 
1215  for(i = 3; i < nEdgeCoeffs; i+=2)
1216  {
1217  signarray[i] = -1;
1218  }
1219  }
1220  }
1221  break;
1222  case 2:
1223  {
1224  for(i = 0; i < nEdgeCoeffs; i++)
1225  {
1226  maparray[i] = nummodes0+i;
1227  }
1228 
1229  if(edgeOrient==eForwards)
1230  {
1231  swap( maparray[0] , maparray[1] );
1232 
1233  for(i = 3; i < nEdgeCoeffs; i+=2)
1234  {
1235  signarray[i] = -1;
1236  }
1237  }
1238  }
1239  break;
1240  case 3:
1241  {
1242  for(i = 0; i < nEdgeCoeffs; i++)
1243  {
1244  maparray[i] = i*nummodes0;
1245  }
1246 
1247  if(edgeOrient==eForwards)
1248  {
1249  swap( maparray[0] , maparray[1] );
1250 
1251  for(i = 3; i < nEdgeCoeffs; i+=2)
1252  {
1253  signarray[i] = -1;
1254  }
1255  }
1256  }
1257  break;
1258  default:
1259  ASSERTL0(false,"eid must be between 0 and 3");
1260  break;
1261  }
1262  }
1263  else if(bType == LibUtilities::eGLL_Lagrange ||
1265  {
1266  switch(eid)
1267  {
1268  case 0:
1269  {
1270  for(i = 0; i < nEdgeCoeffs; i++)
1271  {
1272  maparray[i] = i;
1273  }
1274  }
1275  break;
1276  case 1:
1277  {
1278  for(i = 0; i < nEdgeCoeffs; i++)
1279  {
1280  maparray[i] = (i+1)*nummodes0 - 1;
1281  }
1282  }
1283  break;
1284  case 2:
1285  {
1286  for(i = 0; i < nEdgeCoeffs; i++)
1287  {
1288  maparray[i] = nummodes0*nummodes1 - 1 - i;
1289  }
1290  }
1291  break;
1292  case 3:
1293  {
1294  for(i = 0; i < nEdgeCoeffs; i++)
1295  {
1296  maparray[i] = nummodes0*(nummodes1-1-i);
1297  }
1298  }
1299  break;
1300  default:
1301  ASSERTL0(false,"eid must be between 0 and 3");
1302  break;
1303  }
1304  if(edgeOrient == eBackwards)
1305  {
1306  reverse( maparray.get() , maparray.get()+nEdgeCoeffs );
1307  }
1308  }
1309  else
1310  {
1311  ASSERTL0(false,"Mapping not defined for this type of basis");
1312  }
1313  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
Principle Modified Functions .
Definition: BasisType.h:49
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:54
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:397
Lagrange for SEM basis .
Definition: BasisType.h:53
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdQuadExp::v_GetInteriorMap ( Array< OneD, unsigned int > &  outarray)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 846 of file StdQuadExp.cpp.

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

847  {
848  int i,j;
849  int cnt=0;
850  int nummodes0, nummodes1;
851  int startvalue;
852  if(outarray.num_elements()!=GetNcoeffs()-NumBndryCoeffs())
853  {
854  outarray = Array<OneD, unsigned int>(GetNcoeffs()-NumBndryCoeffs());
855  }
856 
857  nummodes0 = m_base[0]->GetNumModes();
858  nummodes1 = m_base[1]->GetNumModes();
859 
860  const LibUtilities::BasisType Btype0 = GetBasisType(0);
861  const LibUtilities::BasisType Btype1 = GetBasisType(1);
862 
863  switch(Btype1)
864  {
866  startvalue = nummodes0;
867  break;
869  startvalue = 2*nummodes0;
870  break;
871  default:
872  ASSERTL0(0,"Mapping array is not defined for this expansion");
873  break;
874  }
875 
876  switch(Btype0)
877  {
879  startvalue++;
880  break;
882  startvalue+=2;
883  break;
884  default:
885  ASSERTL0(0,"Mapping array is not defined for this expansion");
886  break;
887  }
888 
889  for(i = 0; i < nummodes1-2; i++)
890  {
891  for(j = 0; j < nummodes0-2; j++)
892  {
893  outarray[cnt++]=startvalue+j;
894  }
895  startvalue+=nummodes0;
896  }
897  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
Principle Modified Functions .
Definition: BasisType.h:49
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:131
Lagrange for SEM basis .
Definition: BasisType.h:53
Array< OneD, LibUtilities::BasisSharedPtr > m_base
int Nektar::StdRegions::StdQuadExp::v_GetNedges ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 621 of file StdQuadExp.cpp.

622  {
623  return 4;
624  }
int Nektar::StdRegions::StdQuadExp::v_GetNverts ( ) const
protectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 616 of file StdQuadExp.cpp.

617  {
618  return 4;
619  }
void Nektar::StdRegions::StdQuadExp::v_GetSimplexEquiSpacedConnectivity ( Array< OneD, int > &  conn,
bool  standard = true 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1647 of file StdQuadExp.cpp.

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

1650  {
1651  int np1 = m_base[0]->GetNumPoints();
1652  int np2 = m_base[1]->GetNumPoints();
1653  int np = max(np1,np2);
1654 
1655  conn = Array<OneD, int>(6*(np-1)*(np-1));
1656 
1657  int row = 0;
1658  int rowp1 = 0;
1659  int cnt = 0;
1660  for(int i = 0; i < np-1; ++i)
1661  {
1662  rowp1 += np;
1663  for(int j = 0; j < np-1; ++j)
1664  {
1665  conn[cnt++] = row +j;
1666  conn[cnt++] = row +j+1;
1667  conn[cnt++] = rowp1 +j;
1668 
1669  conn[cnt++] = rowp1 +j+1;
1670  conn[cnt++] = rowp1 +j;
1671  conn[cnt++] = row +j+1;
1672  }
1673  row += np;
1674  }
1675  }
Array< OneD, LibUtilities::BasisSharedPtr > m_base
int Nektar::StdRegions::StdQuadExp::v_GetVertexMap ( int  localVertexId,
bool  useCoeffPacking = false 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 899 of file StdQuadExp.cpp.

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

900  {
901  int localDOF = 0;
902 
903  if(useCoeffPacking == true)
904  {
905  switch(localVertexId)
906  {
907  case 0:
908  {
909  localDOF = 0;
910  }
911  break;
912  case 1:
913  {
915  {
916  localDOF = m_base[0]->GetNumModes()-1;
917  }
918  else
919  {
920  localDOF = 1;
921  }
922  }
923  break;
924  case 2:
925  {
927  {
928  localDOF = m_base[0]->GetNumModes() * (m_base[1]->GetNumModes()-1);
929  }
930  else
931  {
932  localDOF = m_base[0]->GetNumModes();
933  }
934  }
935  break;
936  case 3:
937  {
939  {
940  localDOF = m_base[0]->GetNumModes()*m_base[1]->GetNumModes()-1;
941  }
942  else
943  {
944  localDOF = m_base[0]->GetNumModes()+1;
945  }
946  }
947  break;
948  default:
949  ASSERTL0(false,"eid must be between 0 and 3");
950  break;
951  }
952 
953  }
954  else
955  {
956  switch(localVertexId)
957  {
958  case 0:
959  {
960  localDOF = 0;
961  }
962  break;
963  case 1:
964  {
966  {
967  localDOF = m_base[0]->GetNumModes()-1;
968  }
969  else
970  {
971  localDOF = 1;
972  }
973  }
974  break;
975  case 2:
976  {
978  {
979  localDOF = m_base[0]->GetNumModes()*m_base[1]->GetNumModes()-1;
980  }
981  else
982  {
983  localDOF = m_base[0]->GetNumModes()+1;
984  }
985  }
986  break;
987  case 3:
988  {
990  {
991  localDOF = m_base[0]->GetNumModes() * (m_base[1]->GetNumModes()-1);
992  }
993  else
994  {
995  localDOF = m_base[0]->GetNumModes();
996  }
997  }
998  break;
999  default:
1000  ASSERTL0(false,"eid must be between 0 and 3");
1001  break;
1002  }
1003  }
1004  return localDOF;
1005  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
Lagrange for SEM basis .
Definition: BasisType.h:53
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdQuadExp::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::QuadExp.

Definition at line 1614 of file StdQuadExp.cpp.

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

1617  {
1618  StdQuadExp::v_HelmholtzMatrixOp_MatFree(inarray,outarray,mkey);
1619  }
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
NekDouble Nektar::StdRegions::StdQuadExp::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::QuadExp.

Definition at line 76 of file StdQuadExp.cpp.

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

78  {
79  Array<OneD, const NekDouble> w0 = m_base[0]->GetW();
80  Array<OneD, const NekDouble> w1 = m_base[1]->GetW();
81 
82  return StdExpansion2D::Integral(inarray,w0,w1);
83  }
NekDouble Integral(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &w0, const Array< OneD, const NekDouble > &w1)
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdQuadExp::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*base1 and put into outarray.

$ \begin{array}{rcl} I_{pq} = (\phi_q \phi_q, u) & = & \sum_{i=0}^{nq_0} \sum_{j=0}^{nq_1} \phi_p(\xi_{0,i}) \phi_q(\xi_{1,j}) w^0_i w^1_j u(\xi_{0,i} \xi_{1,j}) \\ & = & \sum_{i=0}^{nq_0} \phi_p(\xi_{0,i}) \sum_{j=0}^{nq_1} \phi_q(\xi_{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_{qi} = \sum_{j=0}^{nq_1} \phi_q(\xi_{1,j}) \tilde{u}_{i,j} = {\bf B_1 U} $ $ I_{pq} = \sum_{i=0}^{nq_0} \phi_p(\xi_{0,i}) f_{qi} = {\bf B_0 F} $

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 379 of file StdQuadExp.cpp.

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

Referenced by v_FwdTrans().

382  {
383  if(m_base[0]->Collocation() && m_base[1]->Collocation())
384  {
385  MultiplyByQuadratureMetric(inarray,outarray);
386  }
387  else
388  {
389  StdQuadExp::v_IProductWRTBase_SumFac(inarray,outarray);
390  }
391  }
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:902
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
Definition: StdQuadExp.cpp:393
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdQuadExp::v_IProductWRTBase_MatOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 422 of file StdQuadExp.cpp.

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

425  {
426  int nq = GetTotPoints();
427  StdMatrixKey iprodmatkey(eIProductWRTBase,DetShapeType(),*this);
428  DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
429 
430  Blas::Dgemv('N',m_ncoeffs,nq,1.0,iprodmat->GetPtr().get(),
431  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
432  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:454
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:684
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:141
void Nektar::StdRegions::StdQuadExp::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::QuadExp.

Definition at line 393 of file StdQuadExp.cpp.

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

Referenced by v_IProductWRTBase().

397  {
398  int nquad0 = m_base[0]->GetNumPoints();
399  int nquad1 = m_base[1]->GetNumPoints();
400  int order0 = m_base[0]->GetNumModes();
401 
402  if(multiplybyweights)
403  {
404  Array<OneD,NekDouble> tmp(nquad0*nquad1+nquad1*order0);
405  Array<OneD,NekDouble> wsp(tmp+nquad0*nquad1);
406 
407  // multiply by integration constants
408  MultiplyByQuadratureMetric(inarray,tmp);
409  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
410  m_base[1]->GetBdata(),
411  tmp,outarray,wsp,true,true);
412  }
413  else
414  {
415  Array<OneD,NekDouble> wsp(nquad1*order0);
416  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
417  m_base[1]->GetBdata(),
418  inarray,outarray,wsp,true,true);
419  }
420  }
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:902
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::StdQuadExp::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 507 of file StdQuadExp.cpp.

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

515  {
516  int nquad0 = m_base[0]->GetNumPoints();
517  int nquad1 = m_base[1]->GetNumPoints();
518  int nmodes0 = m_base[0]->GetNumModes();
519  int nmodes1 = m_base[1]->GetNumModes();
520 
521  bool colldir0 = doCheckCollDir0?(m_base[0]->Collocation()):false;
522  bool colldir1 = doCheckCollDir1?(m_base[1]->Collocation()):false;
523 
524  if(colldir0 && colldir1)
525  {
526  Vmath::Vcopy(m_ncoeffs,inarray.get(),1,outarray.get(),1);
527  }
528  else if(colldir0)
529  {
530  Blas::Dgemm('N', 'N',nmodes0,nmodes1, nquad1,1.0, inarray.get(),
531  nmodes0, base1.get(), nquad1, 0.0,outarray.get(),nmodes0);
532  }
533  else if(colldir1)
534  {
535  Blas::Dgemm('T','N',nmodes0,nquad1,nquad0,1.0,base0.get(),
536  nquad0,inarray.get(),nquad0,0.0,outarray.get(),nmodes0);
537  }
538  else
539  {
540  ASSERTL1(wsp.num_elements()>=nquad1*nmodes0,"Workspace size is not sufficient");
541 
542 #if 1
543  Blas::Dgemm('T','N',nmodes0,nquad1,nquad0,1.0,base0.get(),
544  nquad0,inarray.get(),nquad0,0.0,wsp.get(),nmodes0);
545 
546 #else
547  for(int i = 0; i < nmodes0; ++i)
548  {
549  for(int j = 0; j < nquad1; ++j)
550  {
551  wsp[j*nmodes0+i] = Blas::Ddot(nquad0,
552  base0.get()+i*nquad0,1,
553  inarray.get()+j*nquad0,1);
554  }
555  }
556 #endif
557  Blas::Dgemm('N', 'N',nmodes0,nmodes1, nquad1,1.0, wsp.get(),
558  nmodes0, base1.get(), nquad1, 0.0,outarray.get(),nmodes0);
559  }
560  }
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
#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 Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1038
void Nektar::StdRegions::StdQuadExp::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::QuadExp.

Definition at line 434 of file StdQuadExp.cpp.

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

437  {
438  StdQuadExp::IProductWRTDerivBase_SumFac(dir,inarray,outarray);
439  }
void IProductWRTDerivBase_SumFac(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void Nektar::StdRegions::StdQuadExp::v_IProductWRTDerivBase_MatOp ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 472 of file StdQuadExp.cpp.

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

475  {
476  ASSERTL0((dir==0)||(dir==1),"input dir is out of range");
477 
478  int nq = GetTotPoints();
479  MatrixType mtype;
480 
481  if(dir) // dir == 1
482  {
483  mtype = eIProductWRTDerivBase1;
484  }
485  else // dir == 0
486  {
487  mtype = eIProductWRTDerivBase0;
488  }
489 
490  StdMatrixKey iprodmatkey(mtype,DetShapeType(),*this);
491  DNekMatSharedPtr iprodmat = GetStdMatrix(iprodmatkey);
492 
493  Blas::Dgemv('N',m_ncoeffs,nq,1.0,iprodmat->GetPtr().get(),
494  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
495  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:454
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:684
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:141
void Nektar::StdRegions::StdQuadExp::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::QuadExp.

Definition at line 441 of file StdQuadExp.cpp.

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

444  {
445  ASSERTL0((dir==0)||(dir==1),"input dir is out of range");
446 
447  int nquad0 = m_base[0]->GetNumPoints();
448  int nquad1 = m_base[1]->GetNumPoints();
449  int nqtot = nquad0*nquad1;
450  int order0 = m_base[0]->GetNumModes();
451 
452  Array<OneD,NekDouble> tmp(nqtot+nquad1*order0);
453  Array<OneD,NekDouble> wsp(tmp+nqtot);
454 
455  // multiply by integration constants
456  MultiplyByQuadratureMetric(inarray,tmp);
457 
458  if(dir) // dir == 1
459  {
460  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
461  m_base[1]->GetDbdata(),
462  tmp,outarray,wsp,true,false);
463  }
464  else // dir == 0
465  {
466  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
467  m_base[1]->GetBdata(),
468  tmp,outarray,wsp,false,true);
469  }
470  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:902
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
bool Nektar::StdRegions::StdQuadExp::v_IsBoundaryInteriorExpansion ( )
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 741 of file StdQuadExp.cpp.

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

742  {
743  bool returnval = false;
744 
747  {
750  {
751  returnval = true;
752  }
753  }
754 
755  return returnval;
756  }
Principle Modified Functions .
Definition: BasisType.h:49
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
Lagrange for SEM basis .
Definition: BasisType.h:53
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdQuadExp::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::QuadExp.

Definition at line 1590 of file StdQuadExp.cpp.

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

1594  {
1595  StdQuadExp::v_LaplacianMatrixOp_MatFree(inarray,outarray,mkey);
1596  }
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
void Nektar::StdRegions::StdQuadExp::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::QuadExp.

Definition at line 1598 of file StdQuadExp.cpp.

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

1602  {
1603  StdExpansion::LaplacianMatrixOp_MatFree(k1,k2,inarray,outarray,mkey);
1604  }
void LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void Nektar::StdRegions::StdQuadExp::v_LocCoordToLocCollapsed ( const Array< OneD, const NekDouble > &  xi,
Array< OneD, NekDouble > &  eta 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 567 of file StdQuadExp.cpp.

569  {
570  eta[0] = xi[0];
571  eta[1] = xi[1];
572  }
void Nektar::StdRegions::StdQuadExp::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::QuadExp.

Definition at line 1582 of file StdQuadExp.cpp.

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

1586  {
1587  StdExpansion::MassMatrixOp_MatFree(inarray,outarray,mkey);
1588  }
void MassMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void Nektar::StdRegions::StdQuadExp::v_MultiplyByStdQuadratureMetric ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1622 of file StdQuadExp.cpp.

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

1625  {
1626  int i;
1627  int nquad0 = m_base[0]->GetNumPoints();
1628  int nquad1 = m_base[1]->GetNumPoints();
1629 
1630  const Array<OneD, const NekDouble>& w0 = m_base[0]->GetW();
1631  const Array<OneD, const NekDouble>& w1 = m_base[1]->GetW();
1632 
1633  // multiply by integration constants
1634  for(i = 0; i < nquad1; ++i)
1635  {
1636  Vmath::Vmul(nquad0, inarray.get()+i*nquad0,1,
1637  w0.get(),1,outarray.get()+i*nquad0,1);
1638  }
1639 
1640  for(i = 0; i < nquad0; ++i)
1641  {
1642  Vmath::Vmul(nquad1,outarray.get()+i,nquad0,w1.get(),1,
1643  outarray.get()+i,nquad0);
1644  }
1645  }
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:169
int Nektar::StdRegions::StdQuadExp::v_NumBndryCoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 703 of file StdQuadExp.cpp.

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

704  {
708  "BasisType is not a boundary interior form");
712  "BasisType is not a boundary interior form");
713 
714  return 4 + 2*(GetBasisNumModes(0)-2) + 2*(GetBasisNumModes(1)-2);
715  }
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
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:54
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
Lagrange for SEM basis .
Definition: BasisType.h:53
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
int Nektar::StdRegions::StdQuadExp::v_NumDGBndryCoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 717 of file StdQuadExp.cpp.

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

718  {
722  "BasisType is not a boundary interior form");
726  "BasisType is not a boundary interior form");
727 
728  return 2*GetBasisNumModes(0) + 2*GetBasisNumModes(1);
729  }
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
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:54
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
Lagrange for SEM basis .
Definition: BasisType.h:53
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
void Nektar::StdRegions::StdQuadExp::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.

For quadrilateral region can use the Tensor_Deriv function defined under StdExpansion.

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 95 of file StdQuadExp.cpp.

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

Referenced by v_StdPhysDeriv().

99  {
100  PhysTensorDeriv(inarray, out_d0, out_d1);
101  }
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d0, Array< OneD, NekDouble > &outarray_d1)
Calculate the 2D derivative in the local tensor/collapsed coordinate at the physical points...
void Nektar::StdRegions::StdQuadExp::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::QuadExp.

Definition at line 103 of file StdQuadExp.cpp.

References ASSERTL1, Nektar::NullNekDouble1DArray, and Nektar::StdRegions::StdExpansion2D::PhysTensorDeriv().

106  {
107  switch(dir)
108  {
109  case 0:
110  {
111  PhysTensorDeriv(inarray, outarray, NullNekDouble1DArray);
112  }
113  break;
114  case 1:
115  {
116  PhysTensorDeriv(inarray, NullNekDouble1DArray, outarray);
117  }
118  break;
119  default:
120  {
121  ASSERTL1(false,"input dir is out of range");
122  }
123  break;
124  }
125  }
static Array< OneD, NekDouble > NullNekDouble1DArray
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...
#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::StdQuadExp::v_ReduceOrderCoeffs ( int  numMin,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 1532 of file StdQuadExp.cpp.

References Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eOrtho_A, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::LibUtilities::InterpCoeff2D(), Nektar::StdRegions::StdExpansion::m_base, Vmath::Vcopy(), and Vmath::Zero().

1536  {
1537  int n_coeffs = inarray.num_elements();
1538 
1539 
1540  Array<OneD, NekDouble> coeff(n_coeffs);
1541  Array<OneD, NekDouble> coeff_tmp(n_coeffs,0.0);
1542  Array<OneD, NekDouble> tmp;
1543  Array<OneD, NekDouble> tmp2;
1544 
1545  int nmodes0 = m_base[0]->GetNumModes();
1546  int nmodes1 = m_base[1]->GetNumModes();
1547  int numMax = nmodes0;
1548 
1549  Vmath::Vcopy(n_coeffs,inarray,1,coeff_tmp,1);
1550 
1551  const LibUtilities::PointsKey Pkey0(
1553  const LibUtilities::PointsKey Pkey1(
1555 
1556  LibUtilities::BasisKey b0(m_base[0]->GetBasisType(),nmodes0,Pkey0);
1557  LibUtilities::BasisKey b1(m_base[1]->GetBasisType(),nmodes1,Pkey1);
1558 
1559  LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A,nmodes0,Pkey0);
1560  LibUtilities::BasisKey bortho1(LibUtilities::eOrtho_A,nmodes1,Pkey1);
1561 
1563  b0, b1, coeff_tmp, bortho0, bortho1, coeff);
1564 
1565  Vmath::Zero(n_coeffs,coeff_tmp,1);
1566 
1567  int cnt = 0;
1568  for (int i = 0; i < numMin+1; ++i)
1569  {
1570  Vmath::Vcopy(numMin,
1571  tmp = coeff+cnt,1,
1572  tmp2 = coeff_tmp+cnt,1);
1573 
1574  cnt = i*numMax;
1575  }
1576 
1578  bortho0, bortho1, coeff_tmp, b0, b1, outarray);
1579  }
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
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:359
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void InterpCoeff2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
Definition: InterpCoeff.cpp:73
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1038
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:50
void Nektar::StdRegions::StdQuadExp::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.

Definition at line 127 of file StdQuadExp.cpp.

References v_PhysDeriv().

131  {
132  //PhysTensorDeriv(inarray, out_d0, out_d1);
133  StdQuadExp::v_PhysDeriv(inarray, out_d0, out_d1);
134  }
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: StdQuadExp.cpp:95
void Nektar::StdRegions::StdQuadExp::v_StdPhysDeriv ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 136 of file StdQuadExp.cpp.

References v_PhysDeriv().

139  {
140  //PhysTensorDeriv(inarray, outarray);
141  StdQuadExp::v_PhysDeriv(dir,inarray,outarray);
142 
143  }
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: StdQuadExp.cpp:95
void Nektar::StdRegions::StdQuadExp::v_SVVLaplacianFilter ( Array< OneD, NekDouble > &  array,
const StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::QuadExp.

Definition at line 1477 of file StdQuadExp.cpp.

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

1479  {
1480  // Generate an orthonogal expansion
1481  int qa = m_base[0]->GetNumPoints();
1482  int qb = m_base[1]->GetNumPoints();
1483  int nmodes_a = m_base[0]->GetNumModes();
1484  int nmodes_b = m_base[1]->GetNumModes();
1485  // Declare orthogonal basis.
1486  LibUtilities::PointsKey pa(qa,m_base[0]->GetPointsType());
1487  LibUtilities::PointsKey pb(qb,m_base[1]->GetPointsType());
1488 
1489  LibUtilities::BasisKey Ba(LibUtilities::eOrtho_A,nmodes_a,pa);
1490  LibUtilities::BasisKey Bb(LibUtilities::eOrtho_A,nmodes_b,pb);
1491  StdQuadExp OrthoExp(Ba,Bb);
1492 
1493  Array<OneD, NekDouble> orthocoeffs(OrthoExp.GetNcoeffs());
1494 
1495  //for the "old" implementation
1496  int cutoff = (int) (mkey.GetConstFactor(eFactorSVVCutoffRatio)*min(nmodes_a,nmodes_b));
1497 
1498  //SVV parameters loaded from the .xml case file
1499  NekDouble SvvDiffCoeff = mkey.GetConstFactor(eFactorSVVDiffCoeff);
1500 
1501  // project onto modal space.
1502  OrthoExp.FwdTrans(array,orthocoeffs);
1503 
1504  //counters for scanning through orthocoeffs array
1505  int j, k, cnt = 0;
1506  int nmodes = min(nmodes_a,nmodes_b);
1507 
1508  //------"New" Version August 22nd '13--------------------
1509  for(j = 0; j < nmodes_a; ++j)
1510  {
1511  for(k = 0; k < nmodes_b; ++k)
1512  {
1513  if(j + k >= cutoff) //to filter out only the "high-modes"
1514  {
1515  orthocoeffs[j*nmodes_b+k] *=
1516  (SvvDiffCoeff*exp(-(j+k-nmodes)*(j+k-nmodes)/
1517  ((NekDouble)((j+k-cutoff+1)*
1518  (j+k-cutoff+1)))));
1519  }
1520  else
1521  {
1522  orthocoeffs[j*nmodes_b+k] *= 0.0;
1523  }
1524  cnt++;
1525  }
1526  }
1527 
1528  // backward transform to physical space
1529  OrthoExp.BwdTrans(orthocoeffs,array);
1530  }
Principle Orthogonal Functions .
Definition: BasisType.h:46
double NekDouble
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:216
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdQuadExp::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::QuadExp.

Definition at line 1606 of file StdQuadExp.cpp.

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

1610  {
1611  StdExpansion::WeakDerivMatrixOp_MatFree(i,inarray,outarray,mkey);
1612  }
void WeakDerivMatrixOp_MatFree(const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)