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

Class representing a segment element in reference space. More...

#include <StdSegExp.h>

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

Public Member Functions

 StdSegExp ()
 Default constructor. More...
 
 StdSegExp (const LibUtilities::BasisKey &Ba)
 Constructor using BasisKey class for quadrature points and order definition. More...
 
 StdSegExp (const StdSegExp &T)
 Copy Constructor. More...
 
 ~StdSegExp ()
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion1D
 StdExpansion1D ()
 
 StdExpansion1D (int numcoeffs, const LibUtilities::BasisKey &Ba)
 
 StdExpansion1D (const StdExpansion1D &T)
 
virtual ~StdExpansion1D ()
 
void PhysTensorDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Evaluate the derivative $ d/d{\xi_1} $ at the physical quadrature points given by inarray and return in outarray. More...
 
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion
 StdExpansion ()
 Default Constructor. More...
 
 StdExpansion (const int numcoeffs, const int numbases, const LibUtilities::BasisKey &Ba=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bb=LibUtilities::NullBasisKey, const LibUtilities::BasisKey &Bc=LibUtilities::NullBasisKey)
 Constructor. More...
 
 StdExpansion (const StdExpansion &T)
 Copy Constructor. More...
 
virtual ~StdExpansion ()
 Destructor. More...
 
int GetNumBases () const
 This function returns the number of 1D bases used in the expansion. More...
 
const Array< OneD, const
LibUtilities::BasisSharedPtr > & 
GetBase () const
 This function gets the shared point to basis. More...
 
const
LibUtilities::BasisSharedPtr
GetBasis (int dir) const
 This function gets the shared point to basis in the dir direction. More...
 
int GetNcoeffs (void) const
 This function returns the total number of coefficients used in the expansion. More...
 
int GetTotPoints () const
 This function returns the total number of quadrature points used in the element. More...
 
LibUtilities::BasisType GetBasisType (const int dir) const
 This function returns the type of basis used in the dir direction. More...
 
int GetBasisNumModes (const int dir) const
 This function returns the number of expansion modes in the dir direction. More...
 
int EvalBasisNumModesMax (void) const
 This function returns the maximum number of expansion modes over all local directions. More...
 
LibUtilities::PointsType GetPointsType (const int dir) const
 This function returns the type of quadrature points used in the dir direction. More...
 
int GetNumPoints (const int dir) const
 This function returns the number of quadrature points in the dir direction. More...
 
const Array< OneD, const
NekDouble > & 
GetPoints (const int dir) const
 This function returns a pointer to the array containing the quadrature points in dir direction. More...
 
int GetNverts () const
 This function returns the number of vertices of the expansion domain. More...
 
int GetNedges () const
 This function returns the number of edges of the expansion domain. More...
 
int GetEdgeNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th edge. More...
 
int GetTotalEdgeIntNcoeffs () const
 
int GetEdgeNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th edge. More...
 
int DetCartesianDirOfEdge (const int edge)
 
const LibUtilities::BasisKey DetEdgeBasisKey (const int i) const
 
const LibUtilities::BasisKey DetFaceBasisKey (const int i, const int k) const
 
int GetFaceNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th face. More...
 
int GetFaceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th face. More...
 
int GetFaceIntNcoeffs (const int i) const
 
int GetTotalFaceIntNcoeffs () const
 
int GetTraceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th edge/face. More...
 
LibUtilities::PointsKey GetFacePointsKey (const int i, const int j) const
 
int NumBndryCoeffs (void) const
 
int NumDGBndryCoeffs (void) const
 
LibUtilities::BasisType GetEdgeBasisType (const int i) const
 This function returns the type of expansion basis on the i-th edge. More...
 
const LibUtilities::PointsKey GetNodalPointsKey () const
 This function returns the type of expansion Nodal point type if defined. More...
 
int GetNfaces () const
 This function returns the number of faces of the expansion domain. More...
 
int GetNtrace () const
 Returns the number of trace elements connected to this element. More...
 
LibUtilities::ShapeType DetShapeType () const
 This function returns the shape of the expansion domain. More...
 
boost::shared_ptr< StdExpansionGetStdExp (void) const
 
boost::shared_ptr< StdExpansionGetLinStdExp (void) const
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion ()
 
bool IsNodalNonTensorialExp ()
 
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space. More...
 
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space. More...
 
void FwdTrans_BndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble Integral (const Array< OneD, const NekDouble > &inarray)
 This function integrates the specified function over the domain. More...
 
void FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 This function fills the array outarray with the mode-th mode of the expansion. More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 this function calculates the inner product of a given function f with the different modes of the expansion More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check)
 
void IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
int GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void SetElmtId (const int id)
 Set the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2=NullNekDouble1DArray, Array< OneD, NekDouble > &coords_3=NullNekDouble1DArray)
 this function returns the physical coordinates of the quadrature points of the expansion More...
 
void GetCoord (const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
 given the coordinates of a point of the element in the local collapsed coordinate system, this function calculates the physical coordinates of the point More...
 
DNekMatSharedPtr GetStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr GetStdStaticCondMatrix (const StdMatrixKey &mkey)
 
IndexMapValuesSharedPtr GetIndexMap (const IndexMapKey &ikey)
 
const Array< OneD, const
NekDouble > & 
GetPhysNormals (void)
 
void SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
virtual void SetUpPhysNormals (const int edge)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
void DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
StdRegions::Orientation GetForient (int face)
 
StdRegions::Orientation GetEorient (int edge)
 
StdRegions::Orientation GetPorient (int point)
 
StdRegions::Orientation GetCartesianEorient (int edge)
 
void SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
 
void SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
int CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
NekDouble StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
int GetCoordim ()
 
void GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
void GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
int GetVertexMap (const int localVertexId, bool useCoeffPacking=false)
 
void GetEdgeInteriorMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
 
void GetFaceNumModes (const int fid, const Orientation faceOrient, int &numModes0, int &numModes1)
 
void GetFaceInteriorMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
 
void GetEdgeToElementMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int P=-1)
 
void GetFaceToElementMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int nummodesA=-1, int nummodesB=-1)
 
void GetEdgePhysVals (const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Extract the physical values along edge edge from inarray into outarray following the local edge orientation and point distribution defined by defined in EdgeExp. More...
 
void GetEdgePhysVals (const int edge, const boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetTracePhysVals (const int edge, const boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetVertexPhysVals (const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray)
 
void GetEdgeInterpVals (const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetEdgeQFactors (const int edge, Array< OneD, NekDouble > &outarray)
 Extract the metric factors to compute the contravariant fluxes along edge edge and stores them into outarray following the local edge orientation (i.e. anticlockwise convention). More...
 
void GetFacePhysVals (const int face, const boost::shared_ptr< StdExpansion > &FaceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=eNoOrientation)
 
void GetEdgePhysMap (const int edge, Array< OneD, int > &outarray)
 
void GetFacePhysMap (const int face, Array< OneD, int > &outarray)
 
void MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
DNekMatSharedPtr CreateGeneralMatrix (const StdMatrixKey &mkey)
 this function generates the mass matrix $\mathbf{M}[i][j] = \int \phi_i(\mathbf{x}) \phi_j(\mathbf{x}) d\mathbf{x}$ More...
 
void GeneralMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
DNekMatSharedPtr GenMatrix (const StdMatrixKey &mkey)
 
void PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void PhysDeriv_s (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds)
 
void PhysDeriv_n (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn)
 
void PhysDirectionalDeriv (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &outarray)
 
void StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void AddRobinMassMatrix (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
void AddRobinEdgeContribution (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, Array< OneD, NekDouble > &coeffs)
 
NekDouble PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta. More...
 
const boost::shared_ptr
< SpatialDomains::GeomFactors > & 
GetMetricInfo (void) const
 
virtual int v_GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
virtual const Array< OneD,
const NekDouble > & 
v_GetPhysNormals (void)
 
virtual void v_SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual void v_DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
virtual StdRegions::Orientation v_GetForient (int face)
 
virtual StdRegions::Orientation v_GetEorient (int edge)
 
virtual StdRegions::Orientation v_GetCartesianEorient (int edge)
 
virtual StdRegions::Orientation v_GetPorient (int point)
 
NekDouble Linf (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete $ L_\infty$ error $ |\epsilon|_\infty = \max |u - u_{exact}|$ where $ u_{exact}$ is given by the array sol. More...
 
NekDouble L2 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete $ L_2$ error, $ | \epsilon |_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 dx \right]^{1/2} d\xi_1 $ where $ u_{exact}$ is given by the array sol. More...
 
NekDouble H1 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete $ H^1$ error, $ | \epsilon |^1_{2} = \left [ \int^1_{-1} [u - u_{exact}]^2 + \nabla(u - u_{exact})\cdot\nabla(u - u_{exact})\cdot dx \right]^{1/2} d\xi_1 $ where $ u_{exact}$ is given by the array sol. More...
 
const NormalVectorGetEdgeNormal (const int edge) const
 
void ComputeEdgeNormal (const int edge)
 
void NegateEdgeNormal (const int edge)
 
bool EdgeNormalNegated (const int edge)
 
void ComputeFaceNormal (const int face)
 
void NegateFaceNormal (const int face)
 
bool FaceNormalNegated (const int face)
 
void ComputeVertexNormal (const int vertex)
 
const NormalVectorGetFaceNormal (const int face) const
 
const NormalVectorGetVertexNormal (const int vertex) const
 
const NormalVectorGetSurfaceNormal (const int id) const
 
const LibUtilities::PointsKeyVector GetPointsKeys () const
 
Array< OneD, unsigned int > GetEdgeInverseBoundaryMap (int eid)
 
Array< OneD, unsigned int > GetFaceInverseBoundaryMap (int fid, StdRegions::Orientation faceOrient=eNoOrientation)
 
DNekMatSharedPtr BuildInverseTransformationMatrix (const DNekScalMatSharedPtr &m_transformationmatrix)
 
void PhysInterpToSimplexEquiSpaced (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int npset=-1)
 This function performs an interpolation from the physical space points provided at input into an array of equispaced points which are not the collapsed coordinate. So for a tetrahedron you will only get a tetrahedral number of values. More...
 
void GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true)
 This function provides the connectivity of local simplices (triangles or tets) to connect the equispaced data points provided by PhysInterpToSimplexEquiSpaced. More...
 
void EquiSpacedToCoeffs (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs a projection/interpolation from the equispaced points sometimes used in post-processing onto the coefficient space. More...
 
template<class T >
boost::shared_ptr< T > as ()
 
void IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 

Protected Member Functions

virtual NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray)
 Integrate the physical point list inarray over region and return the value. More...
 
virtual void v_PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 Evaluate the derivative $ d/d{\xi_1} $ at the physical quadrature points given by inarray and return in outarray. 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=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
virtual void v_StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Backward transform from coefficient space given in inarray and evaluate at the physical quadrature points outarray. More...
 
virtual void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Forward transform from physical quadrature space stored in inarray and evaluate the expansion coefficients and store in outarray. More...
 
virtual void v_BwdTrans_SumFac (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_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Inner product of inarray over region with respect to the expansion basis (this)->m_base[0] and return in outarray. More...
 
virtual void v_IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check)
 Inner product of inarray over region with respect to expansion basis base and return in outarray. More...
 
virtual void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
virtual void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &Lcoords, const Array< OneD, const NekDouble > &physvals)
 
virtual void v_LaplacianMatrixOp (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_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 
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 int v_GetNverts () const
 
virtual int v_NumBndryCoeffs () const
 
virtual int v_NumDGBndryCoeffs () const
 
virtual bool v_IsBoundaryInteriorExpansion ()
 
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
virtual LibUtilities::ShapeType v_DetShapeType () const
 Return Shape of region, using ShapeType enum list. i.e. Segment. More...
 
virtual DNekMatSharedPtr v_GenMatrix (const StdMatrixKey &mkey)
 
virtual DNekMatSharedPtr v_CreateStdMatrix (const StdMatrixKey &mkey)
 
virtual void v_GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true)
 
virtual void v_ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion
DNekMatSharedPtr CreateStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr CreateStdStaticCondMatrix (const StdMatrixKey &mkey)
 Create the static condensation of a matrix when using a boundary interior decomposition. More...
 
IndexMapValuesSharedPtr CreateIndexMap (const IndexMapKey &ikey)
 Create an IndexMap which contains mapping information linking any specific element shape with either its boundaries, edges, faces, verteces, etc. More...
 
void BwdTrans_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GeneralMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
 
void LaplacianMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp_MatFree (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void HelmholtzMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
virtual void v_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
 
virtual NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 

Additional Inherited Members

- Protected Attributes inherited from Nektar::StdRegions::StdExpansion1D
std::map< int, NormalVectorm_vertexNormals
 
- Protected Attributes inherited from Nektar::StdRegions::StdExpansion
Array< OneD,
LibUtilities::BasisSharedPtr
m_base
 
int m_elmt_id
 
int m_ncoeffs
 
LibUtilities::NekManager
< StdMatrixKey, DNekMat,
StdMatrixKey::opLess
m_stdMatrixManager
 
LibUtilities::NekManager
< StdMatrixKey, DNekBlkMat,
StdMatrixKey::opLess
m_stdStaticCondMatrixManager
 
LibUtilities::NekManager
< IndexMapKey, IndexMapValues,
IndexMapKey::opLess
m_IndexMapManager
 

Detailed Description

Class representing a segment element in reference space.

All interface of this class sits in StdExpansion class

Definition at line 54 of file StdSegExp.h.

Constructor & Destructor Documentation

Nektar::StdRegions::StdSegExp::StdSegExp ( )

Default constructor.

Definition at line 49 of file StdSegExp.cpp.

50  {
51  }
Nektar::StdRegions::StdSegExp::StdSegExp ( const LibUtilities::BasisKey Ba)

Constructor using BasisKey class for quadrature points and order definition.

Parameters
BaBasisKey class definition containing order and quadrature points.

Definition at line 61 of file StdSegExp.cpp.

61  :
62  StdExpansion(Ba.GetNumModes(), 1, Ba),
63  StdExpansion1D(Ba.GetNumModes(),Ba)
64  {
65  }
StdExpansion()
Default Constructor.
Nektar::StdRegions::StdSegExp::StdSegExp ( const StdSegExp T)

Copy Constructor.

Definition at line 70 of file StdSegExp.cpp.

70  :
71  StdExpansion(T),
73  {
74  }
StdExpansion()
Default Constructor.
Nektar::StdRegions::StdSegExp::~StdSegExp ( )

Definition at line 77 of file StdSegExp.cpp.

78  {
79  }

Member Function Documentation

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

Backward transform from coefficient space given in inarray and evaluate at the physical quadrature points outarray.

Operation can be evaluated as $ u(\xi_{1i}) = \sum_{p=0}^{order-1} \hat{u}_p \phi_p(\xi_{1i}) $ or equivalently $ {\bf u} = {\bf B}^T {\bf \hat{u}} $ where ${\bf B}[i][j] = \phi_i(\xi_{1j}), \mbox{\_coeffs}[p] = {\bf \hat{u}}[p] $

The function takes the coefficient array inarray as input for the transformation

Parameters
inarraythe coeffficients of the expansion
outarraythe resulting array of the values of the function at the physical quadrature points will be stored in the array outarray

Implements Nektar::StdRegions::StdExpansion.

Definition at line 218 of file StdSegExp.cpp.

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

Referenced by v_BwdTrans_SumFac(), v_HelmholtzMatrixOp(), and v_LaplacianMatrixOp().

221  {
222  int nquad = m_base[0]->GetNumPoints();
223 
224  if(m_base[0]->Collocation())
225  {
226  Vmath::Vcopy(nquad, inarray, 1, outarray, 1);
227  }
228  else
229  {
230 
231 #ifdef NEKTAR_USING_DIRECT_BLAS_CALLS
232 
233  Blas::Dgemv('N',nquad,m_base[0]->GetNumModes(),1.0, (m_base[0]->GetBdata()).get(),
234  nquad,&inarray[0],1,0.0,&outarray[0],1);
235 
236 #else //NEKTAR_USING_DIRECT_BLAS_CALLS
237 
238  NekVector<NekDouble> in(m_ncoeffs,inarray,eWrapper);
239  NekVector<NekDouble> out(nquad,outarray,eWrapper);
240  NekMatrix<NekDouble> B(nquad,m_ncoeffs,m_base[0]->GetBdata(),eWrapper);
241  out = B * in;
242 
243 #endif //NEKTAR_USING_DIRECT_BLAS_CALLS
244 
245  }
246  }
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
void Nektar::StdRegions::StdSegExp::v_BwdTrans_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 398 of file StdSegExp.cpp.

References v_BwdTrans().

400  {
401  v_BwdTrans(inarray, outarray);
402  }
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Backward transform from coefficient space given in inarray and evaluate at the physical quadrature po...
Definition: StdSegExp.cpp:218
int Nektar::StdRegions::StdSegExp::v_CalcNumberOfCoefficients ( const std::vector< unsigned int > &  nummodes,
int &  modes_offset 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 622 of file StdSegExp.cpp.

625  {
626  int nmodes = nummodes[modes_offset];
627  modes_offset += 1;
628 
629  return nmodes;
630  }
DNekMatSharedPtr Nektar::StdRegions::StdSegExp::v_CreateStdMatrix ( const StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 702 of file StdSegExp.cpp.

References v_GenMatrix().

703  {
704  return v_GenMatrix(mkey);
705  }
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey)
Definition: StdSegExp.cpp:636
LibUtilities::ShapeType Nektar::StdRegions::StdSegExp::v_DetShapeType ( ) const
protectedvirtual

Return Shape of region, using ShapeType enum list. i.e. Segment.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 84 of file StdSegExp.cpp.

References Nektar::LibUtilities::eSegment.

Referenced by v_FwdTrans(), v_FwdTrans_BndConstrained(), and v_GenMatrix().

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 518 of file StdSegExp.cpp.

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

519  {
520  int nquad = m_base[0]->GetNumPoints();
521  const NekDouble * base = m_base[0]->GetBdata().get();
522 
523  ASSERTL2(mode <= m_ncoeffs,
524  "calling argument mode is larger than total expansion order");
525 
526  Vmath::Vcopy(nquad,(NekDouble *)base+mode*nquad,1, &outarray[0],1);
527  }
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:250
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
void Nektar::StdRegions::StdSegExp::v_FwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Forward transform from physical quadrature space stored in inarray and evaluate the expansion coefficients and store in outarray.

Perform a forward transform using a Galerkin projection by taking the inner product of the physical points and multiplying by the inverse of the mass matrix using the Solve method of the standard matrix container holding the local mass matrix, i.e. $ {\bf \hat{u}} = {\bf M}^{-1} {\bf I} $ where $ {\bf I}[p] = \int^1_{-1} \phi_p(\xi_1) u(\xi_1) d\xi_1 $

This function stores the expansion coefficients calculated by the transformation in the coefficient space array outarray

Parameters
inarrayarray of physical quadrature points to be transformed
outarraythe coeffficients of the expansion

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 302 of file StdSegExp.cpp.

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

Referenced by v_FwdTrans_BndConstrained().

304  {
305  if(m_base[0]->Collocation())
306  {
307  Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
308  }
309  else
310  {
311  v_IProductWRTBase(inarray,outarray);
312 
313  // get Mass matrix inverse
314  StdMatrixKey masskey(eInvMass,v_DetShapeType(),*this);
315  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
316 
317  NekVector<NekDouble> in(m_ncoeffs,outarray,eCopy);
318  NekVector<NekDouble> out(m_ncoeffs,outarray,eWrapper);
319 
320  out = (*matsys)*in;
321  }
322  }
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:706
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Inner product of inarray over region with respect to the expansion basis (this)->m_base[0] and return...
Definition: StdSegExp.cpp:471
virtual LibUtilities::ShapeType v_DetShapeType() const
Return Shape of region, using ShapeType enum list. i.e. Segment.
Definition: StdSegExp.cpp:84
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
void Nektar::StdRegions::StdSegExp::v_FwdTrans_BndConstrained ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 324 of file StdSegExp.cpp.

References ASSERTL0, Nektar::LibUtilities::eGauss_Lagrange, Nektar::LibUtilities::eGLL_Lagrange, Nektar::StdRegions::eMass, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::GetVertexMap(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::StdRegions::StdExpansion::m_stdStaticCondMatrixManager, Nektar::StdRegions::StdExpansion::MassMatrixOp(), v_DetShapeType(), v_FwdTrans(), v_IProductWRTBase(), Vmath::Vcopy(), and Vmath::Vsub().

327  {
328  if(m_base[0]->Collocation())
329  {
330  Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
331  }
332  else
333  {
334  int nInteriorDofs = m_ncoeffs-2;
335  int offset;
336 
337  switch(m_base[0]->GetBasisType())
338  {
340  {
341  offset = 1;
342  }
343  break;
345  {
346  nInteriorDofs = m_ncoeffs;
347  offset = 0;
348  }
349  break;
352  {
353  offset = 2;
354  }
355  break;
356  default:
357  ASSERTL0(false,"This type of FwdTrans is not defined for this expansion type");
358  }
359 
360  fill(outarray.get(), outarray.get()+m_ncoeffs, 0.0 );
361 
363  {
364  outarray[GetVertexMap(0)] = inarray[0];
365  outarray[GetVertexMap(1)] = inarray[m_base[0]->GetNumPoints()-1];
366 
367  if(m_ncoeffs>2)
368  {
369  // ideally, we would like to have tmp0 to be replaced by
370  // outarray (currently MassMatrixOp does not allow aliasing)
371  Array<OneD, NekDouble> tmp0(m_ncoeffs);
372  Array<OneD, NekDouble> tmp1(m_ncoeffs);
373 
374  StdMatrixKey masskey(eMass,v_DetShapeType(),*this);
375  MassMatrixOp(outarray,tmp0,masskey);
376  v_IProductWRTBase(inarray,tmp1);
377 
378  Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
379 
380  // get Mass matrix inverse (only of interior DOF)
381  DNekMatSharedPtr matsys =
382  (m_stdStaticCondMatrixManager[masskey])-> GetBlock(1,1);
383 
384  Blas::Dgemv('N',nInteriorDofs,nInteriorDofs,1.0,
385  &(matsys->GetPtr())[0],nInteriorDofs,tmp1.get()
386  +offset,1,0.0,outarray.get()+offset,1);
387  }
388  }
389  else
390  {
391  StdSegExp::v_FwdTrans(inarray, outarray);
392  }
393  }
394 
395  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:976
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Forward transform from physical quadrature space stored in inarray and evaluate the expansion coeffic...
Definition: StdSegExp.cpp:302
Principle Modified Functions .
Definition: BasisType.h:49
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:54
LibUtilities::NekManager< StdMatrixKey, DNekBlkMat, StdMatrixKey::opLess > m_stdStaticCondMatrixManager
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
int GetVertexMap(const int localVertexId, bool useCoeffPacking=false)
Definition: StdExpansion.h:824
Principle Modified Functions .
Definition: BasisType.h:50
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Inner product of inarray over region with respect to the expansion basis (this)->m_base[0] and return...
Definition: StdSegExp.cpp:471
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:343
virtual LibUtilities::ShapeType v_DetShapeType() const
Return Shape of region, using ShapeType enum list. i.e. Segment.
Definition: StdSegExp.cpp:84
Lagrange for SEM basis .
Definition: BasisType.h:53
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
DNekMatSharedPtr Nektar::StdRegions::StdSegExp::v_GenMatrix ( const StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 636 of file StdSegExp.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::StdRegions::StdMatrixKey::ConstFactorExists(), Nektar::StdRegions::StdExpansion::CreateGeneralMatrix(), Nektar::StdRegions::eFactorConst, Nektar::LibUtilities::eFourier, Nektar::StdRegions::eFwdTrans, Nektar::StdRegions::eInvMass, Nektar::StdRegions::eIProductWRTBase, Nektar::StdRegions::eMass, Nektar::StdRegions::ePhysInterpToEquiSpaced, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdMatrixKey::GetMatrixType(), Nektar::LibUtilities::StdSegData::getNumberOfCoefficients(), Nektar::StdRegions::StdExpansion::GetStdMatrix(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, v_DetShapeType(), and Vmath::Vcopy().

Referenced by v_CreateStdMatrix().

637  {
638  DNekMatSharedPtr Mat;
639  MatrixType mattype;
640 
641  switch(mattype = mkey.GetMatrixType())
642  {
644  {
645  int nq = m_base[0]->GetNumPoints();
646 
647  // take definition from key
648  if(mkey.ConstFactorExists(eFactorConst))
649  {
650  nq = (int) mkey.GetConstFactor(eFactorConst);
651  }
652 
655  Array<OneD, NekDouble > coords (1);
656  DNekMatSharedPtr I ;
658  AllocateSharedPtr(neq, nq);
659 
660  for(int i = 0; i < neq; ++i)
661  {
662  coords[0] = -1.0 + 2*i/(NekDouble)(neq-1);
663  I = m_base[0]->GetI(coords);
664  Vmath::Vcopy(nq, I->GetRawPtr(), 1,
665  Mat->GetRawPtr()+i,neq);
666  }
667  }
668  break;
669  case eFwdTrans:
670  {
672  StdMatrixKey iprodkey(eIProductWRTBase,v_DetShapeType(),*this);
673  DNekMat &Iprod = *GetStdMatrix(iprodkey);
674  StdMatrixKey imasskey(eInvMass,v_DetShapeType(),*this);
675  DNekMat &Imass = *GetStdMatrix(imasskey);
676 
677  (*Mat) = Imass*Iprod;
678  }
679  break;
680  default:
681  {
683 
684  if(mattype == eMass)
685  {
686  // For Fourier basis set the imaginary component
687  // of mean mode to have a unit diagonal component
688  // in mass matrix
690  {
691  (*Mat)(1,1) = 1.0;
692  }
693  }
694  }
695  break;
696  }
697 
698  return Mat;
699  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
Fourier Expansion .
Definition: BasisType.h:52
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:706
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
virtual LibUtilities::ShapeType v_DetShapeType() const
Return Shape of region, using ShapeType enum list. i.e. Segment.
Definition: StdSegExp.cpp:84
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
void Nektar::StdRegions::StdSegExp::v_GetBoundaryMap ( Array< OneD, unsigned int > &  outarray)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 712 of file StdSegExp.cpp.

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

713  {
714  if(outarray.num_elements() != NumBndryCoeffs())
715  {
716  outarray = Array<OneD, unsigned int>(NumBndryCoeffs());
717  }
718  const LibUtilities::BasisType Btype = GetBasisType(0);
719  int nummodes = m_base[0]->GetNumModes();
720 
721  outarray[0] = 0;
722 
723  switch(Btype)
724  {
729  outarray[1]= nummodes-1;
730  break;
733  outarray[1] = 1;
734  break;
735  default:
736  ASSERTL0(0,"Mapping array is not defined for this expansion");
737  break;
738  }
739  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
Principle Modified Functions .
Definition: BasisType.h:49
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:54
Fourier Expansion .
Definition: BasisType.h:52
Chebyshev Polynomials .
Definition: BasisType.h:56
Principle Modified Functions .
Definition: BasisType.h:50
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
Lagrange for SEM basis .
Definition: BasisType.h:53
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdSegExp::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::SegExp.

Definition at line 590 of file StdSegExp.cpp.

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

594  {
595  Blas::Dcopy(GetNumPoints(0),(m_base[0]->GetZ()).get(),
596  1,&coords_0[0],1);
597  }
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::StdSegExp::v_GetInteriorMap ( Array< OneD, unsigned int > &  outarray)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 741 of file StdSegExp.cpp.

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

742  {
743  int i;
744  if(outarray.num_elements()!=GetNcoeffs()-NumBndryCoeffs())
745  {
746  outarray = Array<OneD, unsigned int>(GetNcoeffs()-NumBndryCoeffs());
747  }
748  const LibUtilities::BasisType Btype = GetBasisType(0);
749 
750  switch(Btype)
751  {
756  for(i = 0 ; i < GetNcoeffs()-2;i++)
757  {
758  outarray[i] = i+1;
759  }
760  break;
763  for(i = 0 ; i < GetNcoeffs()-2;i++)
764  {
765  outarray[i] = i+2;
766  }
767  break;
768  default:
769  ASSERTL0(0,"Mapping array is not defined for this expansion");
770  break;
771  }
772  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
Principle Modified Functions .
Definition: BasisType.h:49
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:54
Fourier Expansion .
Definition: BasisType.h:52
Chebyshev Polynomials .
Definition: BasisType.h:56
Principle Modified Functions .
Definition: BasisType.h:50
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:131
Lagrange for SEM basis .
Definition: BasisType.h:53
int Nektar::StdRegions::StdSegExp::v_GetNverts ( ) const
protectedvirtual

Implements Nektar::StdRegions::StdExpansion.

Definition at line 607 of file StdSegExp.cpp.

608  {
609  return 2;
610  }
void Nektar::StdRegions::StdSegExp::v_GetSimplexEquiSpacedConnectivity ( Array< OneD, int > &  conn,
bool  standard = true 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 789 of file StdSegExp.cpp.

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

792  {
793  int np = m_base[0]->GetNumPoints();
794 
795  conn = Array<OneD, int>(2*(np-1));
796  int cnt = 0;
797  for(int i = 0; i < np-1; ++i)
798  {
799  conn[cnt++] = i;
800  conn[cnt++] = i+1;
801  }
802  }
Array< OneD, LibUtilities::BasisSharedPtr > m_base
int Nektar::StdRegions::StdSegExp::v_GetVertexMap ( int  localVertexId,
bool  useCoeffPacking = false 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 774 of file StdSegExp.cpp.

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

775  {
776  ASSERTL0((localVertexId==0)||(localVertexId==1),"local vertex id"
777  "must be between 0 or 1");
778 
779  int localDOF = localVertexId;
780 
782  (localVertexId==1) )
783  {
784  localDOF = m_base[0]->GetNumModes()-1;
785  }
786  return localDOF;
787  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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::StdSegExp::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::SegExp.

Definition at line 555 of file StdSegExp.cpp.

References Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, v_BwdTrans(), v_IProductWRTBase(), and v_PhysDeriv().

559  {
560  int nquad = m_base[0]->GetNumPoints();
561 
562  Array<OneD,NekDouble> physValues(nquad);
563  Array<OneD,NekDouble> dPhysValuesdx(nquad);
564  Array<OneD,NekDouble> wsp(m_ncoeffs);
565 
566  v_BwdTrans(inarray,physValues);
567 
568  // mass matrix operation
569  v_IProductWRTBase((m_base[0]->GetBdata()),physValues,wsp,1);
570 
571  // Laplacian matrix operation
572  v_PhysDeriv(physValues,dPhysValuesdx);
573  v_IProductWRTBase(m_base[0]->GetDbdata(),dPhysValuesdx,outarray,1);
574  Blas::Daxpy(m_ncoeffs, mkey.GetConstFactor(eFactorLambda), wsp.get(), 1, outarray.get(), 1);
575  }
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Backward transform from coefficient space given in inarray and evaluate at the physical quadrature po...
Definition: StdSegExp.cpp:218
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Inner product of inarray over region with respect to the expansion basis (this)->m_base[0] and return...
Definition: StdSegExp.cpp:471
Array< OneD, LibUtilities::BasisSharedPtr > m_base
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Evaluate the derivative at the physical quadrature points given by inarray and return in outarray...
Definition: StdSegExp.cpp:155
NekDouble Nektar::StdRegions::StdSegExp::v_Integral ( const Array< OneD, const NekDouble > &  inarray)
protectedvirtual

Integrate the physical point list inarray over region and return the value.

Parameters
inarraydefinition of function to be integrated evauluated at quadrature point of expansion.
Returns
returns $\int^1_{-1} u(\xi_1)d \xi_1 $ where $inarray[i] = u(\xi_{1i}) $

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 122 of file StdSegExp.cpp.

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

123  {
124  NekDouble Int = 0.0;
125  int nquad0 = m_base[0]->GetNumPoints();
126  Array<OneD, NekDouble> tmp(nquad0);
127  Array<OneD, const NekDouble> z = m_base[0]->GetZ();
128  Array<OneD, const NekDouble> w0 = m_base[0]->GetW();
129 
130  // multiply by integration constants
131  Vmath::Vmul(nquad0, inarray, 1, w0, 1, tmp, 1);
132 
133  Int = Vmath::Vsum(nquad0, tmp, 1);
134 
135  return Int;
136  }
double NekDouble
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:737
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:183
void Nektar::StdRegions::StdSegExp::v_IProductWRTBase ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Inner product of inarray over region with respect to the expansion basis (this)->m_base[0] and return in outarray.

Wrapper call to StdSegExp::IProductWRTBase

Parameters
inarrayarray of function values evaluated at the physical collocation points
outarraythe values of the inner product with respect to each basis over region will be stored in the array outarray as output of the function

Implements Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 471 of file StdSegExp.cpp.

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

Referenced by v_FwdTrans(), v_FwdTrans_BndConstrained(), v_HelmholtzMatrixOp(), v_IProductWRTDerivBase(), and v_LaplacianMatrixOp().

473  {
474  v_IProductWRTBase(m_base[0]->GetBdata(),inarray,outarray,1);
475  }
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Inner product of inarray over region with respect to the expansion basis (this)->m_base[0] and return...
Definition: StdSegExp.cpp:471
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::StdRegions::StdSegExp::v_IProductWRTBase ( const Array< OneD, const NekDouble > &  base,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
int  coll_check 
)
protectedvirtual

Inner product of inarray over region with respect to expansion basis base and return in outarray.

Calculate $ I[p] = \int^{1}_{-1} \phi_p(\xi_1) u(\xi_1) d\xi_1 = \sum_{i=0}^{nq-1} \phi_p(\xi_{1i}) u(\xi_{1i}) w_i $ where $ outarray[p] = I[p], inarray[i] = u(\xi_{1i}), base[p*nq+i] = \phi_p(\xi_{1i}) $.

Parameters
basean array defining the local basis for the inner product usually passed from Basis->GetBdata() or Basis->GetDbdata()
inarrayphysical point array of function to be integrated $ u(\xi_1) $
coll_checkflag to identify when a Basis->Collocation() call should be performed to see if this is a GLL_Lagrange basis with a collocation property. (should be set to 0 if taking the inner product with respect to the derivative of basis)
outarraythe values of the inner product with respect to each basis over region will be stored in the array outarray as output of the function

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 433 of file StdSegExp.cpp.

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

438  {
439  int nquad = m_base[0]->GetNumPoints();
440  Array<OneD, NekDouble> tmp(nquad);
441  Array<OneD, const NekDouble> w = m_base[0]->GetW();
442 
443  Vmath::Vmul(nquad, inarray, 1, w, 1, tmp, 1);
444 
445  /* Comment below was a bug for collocated basis
446  if(coll_check&&m_base[0]->Collocation())
447  {
448  Vmath::Vcopy(nquad, tmp, 1, outarray, 1);
449  }
450  else
451  {
452  Blas::Dgemv('T',nquad,m_ncoeffs,1.0,base.get(),nquad,
453  &tmp[0],1,0.0,outarray.get(),1);
454  }*/
455 
456  // Correct implementation
457  Blas::Dgemv('T',nquad,m_ncoeffs,1.0,base.get(),nquad,
458  &tmp[0],1,0.0,outarray.get(),1);
459  }
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:183
void Nektar::StdRegions::StdSegExp::v_IProductWRTBase_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
bool  multiplybyweights = true 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 486 of file StdSegExp.cpp.

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

490  {
491  int nquad = m_base[0]->GetNumPoints();
492  Array<OneD, NekDouble> tmp(nquad);
493  Array<OneD, const NekDouble> w = m_base[0]->GetW();
494  Array<OneD, const NekDouble> base = m_base[0]->GetBdata();
495 
496  if(multiplybyweights)
497  {
498  Vmath::Vmul(nquad, inarray, 1, w, 1, tmp, 1);
499 
500  Blas::Dgemv('T',nquad,m_ncoeffs,1.0,base.get(),nquad,
501  &tmp[0],1,0.0,outarray.get(),1);
502  }
503  else
504  {
505  Blas::Dgemv('T',nquad,m_ncoeffs,1.0,base.get(),nquad,
506  &inarray[0],1,0.0,outarray.get(),1);
507  }
508  }
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:183
void Nektar::StdRegions::StdSegExp::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::SegExp.

Definition at line 477 of file StdSegExp.cpp.

References ASSERTL1, Nektar::StdRegions::StdExpansion::m_base, and v_IProductWRTBase().

481  {
482  ASSERTL1(dir >= 0 && dir < 1,"input dir is out of range");
483  v_IProductWRTBase(m_base[0]->GetDbdata(),inarray,outarray,1);
484  }
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Inner product of inarray over region with respect to the expansion basis (this)->m_base[0] and return...
Definition: StdSegExp.cpp:471
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
Array< OneD, LibUtilities::BasisSharedPtr > m_base
bool Nektar::StdRegions::StdSegExp::v_IsBoundaryInteriorExpansion ( )
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 89 of file StdSegExp.cpp.

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

90  {
91 
92  bool returnval = false;
93 
95  {
96  returnval = true;
97  }
98 
100  {
101  returnval = true;
102  }
103 
104  return returnval;
105  }
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::StdSegExp::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::SegExp.

Definition at line 537 of file StdSegExp.cpp.

References Nektar::StdRegions::StdExpansion::m_base, v_BwdTrans(), v_IProductWRTBase(), and v_PhysDeriv().

541  {
542  int nquad = m_base[0]->GetNumPoints();
543 
544  Array<OneD,NekDouble> physValues(nquad);
545  Array<OneD,NekDouble> dPhysValuesdx(nquad);
546 
547  v_BwdTrans(inarray,physValues);
548 
549  // Laplacian matrix operation
550  v_PhysDeriv(physValues,dPhysValuesdx);
551  v_IProductWRTBase(m_base[0]->GetDbdata(),dPhysValuesdx,outarray,1);
552  }
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Backward transform from coefficient space given in inarray and evaluate at the physical quadrature po...
Definition: StdSegExp.cpp:218
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Inner product of inarray over region with respect to the expansion basis (this)->m_base[0] and return...
Definition: StdSegExp.cpp:471
Array< OneD, LibUtilities::BasisSharedPtr > m_base
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
Evaluate the derivative at the physical quadrature points given by inarray and return in outarray...
Definition: StdSegExp.cpp:155
void Nektar::StdRegions::StdSegExp::v_MultiplyByStdQuadratureMetric ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 578 of file StdSegExp.cpp.

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

581  {
582  int nquad0 = m_base[0]->GetNumPoints();
583 
584  const Array<OneD, const NekDouble>& w0 = m_base[0]->GetW();
585 
586  Vmath::Vmul(nquad0, inarray.get(),1,
587  w0.get(),1,outarray.get(),1);
588  }
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:183
int Nektar::StdRegions::StdSegExp::v_NumBndryCoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 612 of file StdSegExp.cpp.

613  {
614  return 2;
615  }
int Nektar::StdRegions::StdSegExp::v_NumDGBndryCoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 617 of file StdSegExp.cpp.

618  {
619  return 2;
620  }
void Nektar::StdRegions::StdSegExp::v_PhysDeriv ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d0,
Array< OneD, NekDouble > &  out_d1 = NullNekDouble1DArray,
Array< OneD, NekDouble > &  out_d2 = NullNekDouble1DArray 
)
protectedvirtual

Evaluate the derivative $ d/d{\xi_1} $ at the physical quadrature points given by inarray and return in outarray.

This is a wrapper around StdExpansion1D::Tensor_Deriv

Parameters
inarrayarray of a function evaluated at the quadrature points
outarraythe resulting array of the derivative $ du/d_{\xi_1}|_{\xi_{1i}} $ will be stored in the array outarra

Reimplemented from Nektar::StdRegions::StdExpansion.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 155 of file StdSegExp.cpp.

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

Referenced by v_HelmholtzMatrixOp(), and v_LaplacianMatrixOp().

159  {
160  PhysTensorDeriv(inarray,out_d0);
161  }
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Evaluate the derivative at the physical quadrature points given by inarray and return in outarray...
void Nektar::StdRegions::StdSegExp::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::SegExp.

Definition at line 164 of file StdSegExp.cpp.

References ASSERTL1, and Nektar::StdRegions::StdExpansion1D::PhysTensorDeriv().

167  {
168  ASSERTL1(dir==0,"input dir is out of range");
169  PhysTensorDeriv(inarray,outarray);
170  // PhysDeriv(inarray, outarray);
171  }
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Evaluate the derivative at the physical quadrature points given by inarray and return in outarray...
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
NekDouble Nektar::StdRegions::StdSegExp::v_PhysEvaluate ( const Array< OneD, const NekDouble > &  Lcoords,
const Array< OneD, const NekDouble > &  physvals 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion1D.

Reimplemented in Nektar::LocalRegions::SegExp.

Definition at line 530 of file StdSegExp.cpp.

References Nektar::StdRegions::StdExpansion1D::v_PhysEvaluate().

533  {
534  return StdExpansion1D::v_PhysEvaluate(coords, physvals);
535  }
virtual NekDouble v_PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
void Nektar::StdRegions::StdSegExp::v_ReduceOrderCoeffs ( int  numMin,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 248 of file StdSegExp.cpp.

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

252  {
253  int n_coeffs = inarray.num_elements();
254 
255  Array<OneD, NekDouble> coeff(n_coeffs);
256  Array<OneD, NekDouble> coeff_tmp(n_coeffs,0.0);
257  Array<OneD, NekDouble> tmp;
258  Array<OneD, NekDouble> tmp2;
259 
260  int nmodes0 = m_base[0]->GetNumModes();
261 
262  Vmath::Vcopy(n_coeffs,inarray,1,coeff_tmp,1);
263 
264  const LibUtilities::PointsKey Pkey0(
266 
267  LibUtilities::BasisKey b0(m_base[0]->GetBasisType(),nmodes0,Pkey0);
268 
269  LibUtilities::BasisKey bortho0(LibUtilities::eOrtho_A,nmodes0,Pkey0);
270 
272  b0, coeff_tmp, bortho0, coeff);
273 
274  Vmath::Zero(n_coeffs,coeff_tmp,1);
275 
276  Vmath::Vcopy(numMin,
277  tmp = coeff,1,
278  tmp2 = coeff_tmp,1);
279 
281  bortho0, coeff_tmp, b0, outarray);
282  }
void InterpCoeff1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
Definition: InterpCoeff.cpp:47
Principle Orthogonal Functions .
Definition: BasisType.h:46
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:373
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:52
void Nektar::StdRegions::StdSegExp::v_StdPhysDeriv ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d0,
Array< OneD, NekDouble > &  out_d1 = NullNekDouble1DArray,
Array< OneD, NekDouble > &  out_d2 = NullNekDouble1DArray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 173 of file StdSegExp.cpp.

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

178  {
179  PhysTensorDeriv(inarray,out_d0);
180  // PhysDeriv(inarray, out_d0);
181  }
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Evaluate the derivative at the physical quadrature points given by inarray and return in outarray...
void Nektar::StdRegions::StdSegExp::v_StdPhysDeriv ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 183 of file StdSegExp.cpp.

References ASSERTL1, and Nektar::StdRegions::StdExpansion1D::PhysTensorDeriv().

187  {
188  ASSERTL1(dir==0,"input dir is out of range");
189  PhysTensorDeriv(inarray,outarray);
190  // PhysDeriv(inarray, outarray);
191  }
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Evaluate the derivative at the physical quadrature points given by inarray and return in outarray...
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228