Nektar++
Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes | List of all members
Nektar::LocalRegions::SegExp Class Reference

#include <SegExp.h>

Inheritance diagram for Nektar::LocalRegions::SegExp:
Inheritance graph
[legend]
Collaboration diagram for Nektar::LocalRegions::SegExp:
Collaboration graph
[legend]

Public Member Functions

 SegExp (const LibUtilities::BasisKey &Ba, const SpatialDomains::Geometry1DSharedPtr &geom)
 Constructor using BasisKey class for quadrature points and order definition. More...
 
 SegExp (const SegExp &S)
 Copy Constructor. More...
 
 ~SegExp ()
 
- Public Member Functions inherited from Nektar::StdRegions::StdSegExp
 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::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 void v_SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
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 StdRegions::Orientation v_GetForient (int face)
 
virtual StdRegions::Orientation v_GetEorient (int edge)
 
virtual StdRegions::Orientation v_GetCartesianEorient (int edge)
 
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)
 
- Public Member Functions inherited from Nektar::LocalRegions::Expansion1D
 Expansion1D (SpatialDomains::Geometry1DSharedPtr pGeom)
 
virtual ~Expansion1D ()
 
void AddNormTraceInt (const int dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
Expansion2DSharedPtr GetLeftAdjacentElementExp () const
 
Expansion2DSharedPtr GetRightAdjacentElementExp () const
 
int GetLeftAdjacentElementEdge () const
 
int GetRightAdjacentElementEdge () const
 
void SetAdjacentElementExp (int edge, Expansion2DSharedPtr &e)
 
void AddHDGHelmholtzTraceTerms (const NekDouble tau, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
SpatialDomains::Geometry1DSharedPtr GetGeom1D () const
 
- Public Member Functions inherited from Nektar::LocalRegions::Expansion
 Expansion (SpatialDomains::GeometrySharedPtr pGeom)
 
 Expansion (const Expansion &pSrc)
 
virtual ~Expansion ()
 
DNekScalMatSharedPtr GetLocMatrix (const LocalRegions::MatrixKey &mkey)
 
DNekScalMatSharedPtr GetLocMatrix (const StdRegions::MatrixType mtype, const StdRegions::ConstFactorMap &factors=StdRegions::NullConstFactorMap, const StdRegions::VarCoeffMap &varcoeffs=StdRegions::NullVarCoeffMap)
 
SpatialDomains::GeometrySharedPtr GetGeom () const
 
void Reset ()
 
virtual const SpatialDomains::GeomFactorsSharedPtrv_GetMetricInfo () const
 
DNekMatSharedPtr BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
 
DNekMatSharedPtr BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd)
 
void AddEdgeNormBoundaryInt (const int edge, const boost::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
void AddEdgeNormBoundaryInt (const int edge, const boost::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
void AddFaceNormBoundaryInt (const int face, const boost::shared_ptr< Expansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
void DGDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &coeffs, Array< OneD, NekDouble > &outarray)
 

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_PhysDeriv_s (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds)
 Evaluate the derivative along a line: $ d/ds=\frac{spacedim}{||tangent||}d/d{\xi} $. The derivative is calculated performing the product $ du/d{s}=\nabla u \cdot tangent $. More...
 
virtual void v_PhysDeriv_n (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn)
 Evaluate the derivative normal to a line: $ d/dn=\frac{spacedim}{||normal||}d/d{\xi} $. The derivative is calculated performing the product $ du/d{s}=\nabla u \cdot normal $. 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_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)->_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_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, 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 NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
virtual NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &coord, const Array< OneD, const NekDouble > &physvals)
 
virtual void v_GetCoord (const Array< OneD, const NekDouble > &Lcoords, Array< OneD, NekDouble > &coords)
 
virtual void v_GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3)
 
virtual void v_GetVertexPhysVals (const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray)
 
virtual StdRegions::StdExpansionSharedPtr v_GetStdExp (void) const
 
virtual int v_GetCoordim ()
 
virtual void v_SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
 
virtual void v_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual int v_GetNumPoints (const int dir) const
 
virtual int v_GetNcoeffs (void) const
 
virtual const LibUtilities::BasisSharedPtrv_GetBasis (int dir) const
 
virtual int v_NumBndryCoeffs () const
 
virtual int v_NumDGBndryCoeffs () const
 
virtual void v_ComputeVertexNormal (const int vertex)
 
virtual StdRegions::Orientation v_GetPorient (int point)
 
virtual SpatialDomains::GeomType v_MetricInfoType ()
 
virtual void v_ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int mode_offset, NekDouble *coeffs)
 Unpack data from input file assuming it comes from. More...
 
virtual const Array< OneD, const NekDouble > & v_GetPhysNormals (void)
 
virtual void v_LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
 
virtual void v_HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
 
virtual DNekMatSharedPtr v_GenMatrix (const StdRegions::StdMatrixKey &mkey)
 
DNekScalMatSharedPtr CreateMatrix (const MatrixKey &mkey)
 
virtual DNekMatSharedPtr v_CreateStdMatrix (const StdRegions::StdMatrixKey &mkey)
 
DNekScalBlkMatSharedPtr CreateStaticCondMatrix (const MatrixKey &mkey)
 
virtual DNekScalMatSharedPtr v_GetLocMatrix (const MatrixKey &mkey)
 
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const MatrixKey &mkey)
 
void v_DropLocStaticCondMatrix (const MatrixKey &mkey)
 
- Protected Member Functions inherited from Nektar::StdRegions::StdSegExp
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_BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
virtual void v_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_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 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...
 
- 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)
 
- Protected Member Functions inherited from Nektar::LocalRegions::Expansion1D
virtual void v_AddRobinMassMatrix (const int vert, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
virtual void v_AddRobinEdgeContribution (const int vert, const Array< OneD, const NekDouble > &primCoeffs, Array< OneD, NekDouble > &coeffs)
 
- Protected Member Functions inherited from Nektar::LocalRegions::Expansion
void ComputeLaplacianMetric ()
 
void ComputeQuadratureMetric ()
 
virtual void v_MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_ComputeLaplacianMetric ()
 
virtual DNekMatSharedPtr v_BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
 
virtual DNekMatSharedPtr v_BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd)
 
virtual void v_AddEdgeNormBoundaryInt (const int edge, const boost::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddEdgeNormBoundaryInt (const int edge, const boost::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddFaceNormBoundaryInt (const int face, const boost::shared_ptr< Expansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
virtual void v_DGDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &coeffs, Array< OneD, NekDouble > &outarray)
 

Private Member Functions

 SegExp ()
 
void ReverseCoeffsAndSign (const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Reverse the coefficients in a boundary interior expansion this routine is of use when we need the segment coefficients corresponding to a expansion in the reverse coordinate direction. More...
 
void MultiplyByElmtInvMass (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 

Private Attributes

LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLessm_matrixManager
 
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLessm_staticCondMatrixManager
 

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::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
 
- Protected Attributes inherited from Nektar::LocalRegions::Expansion
SpatialDomains::GeometrySharedPtr m_geom
 
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
 
MetricMap m_metrics
 

Detailed Description

Defines a Segment local expansion.

Definition at line 52 of file SegExp.h.

Constructor & Destructor Documentation

Nektar::LocalRegions::SegExp::SegExp ( const LibUtilities::BasisKey Ba,
const SpatialDomains::Geometry1DSharedPtr geom 
)

Constructor using BasisKey class for quadrature points and order definition.

Parameters
BaBasis key of segment expansion.
geomDescription of geometry.

Definition at line 57 of file SegExp.cpp.

58  :
59  StdExpansion(Ba.GetNumModes(), 1, Ba),
60  StdExpansion1D(Ba.GetNumModes(), Ba),
61  StdRegions::StdSegExp(Ba),
62  Expansion(geom),
63  Expansion1D(geom),
65  boost::bind(&SegExp::CreateMatrix, this, _1),
66  std::string("SegExpMatrix")),
68  boost::bind(&SegExp::CreateStaticCondMatrix, this, _1),
69  std::string("SegExpStaticCondMatrix"))
70  {
71  }
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: SegExp.h:237
Expansion(SpatialDomains::GeometrySharedPtr pGeom)
Definition: Expansion.cpp:46
Expansion1D(SpatialDomains::Geometry1DSharedPtr pGeom)
Definition: Expansion1D.h:63
DNekScalBlkMatSharedPtr CreateStaticCondMatrix(const MatrixKey &mkey)
Definition: SegExp.cpp:1411
DNekScalMatSharedPtr CreateMatrix(const MatrixKey &mkey)
Definition: SegExp.cpp:1177
StdExpansion()
Default Constructor.
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: SegExp.h:235
Nektar::LocalRegions::SegExp::SegExp ( const SegExp S)

Copy Constructor.

Parameters
SExisting segment to duplicate.

Definition at line 78 of file SegExp.cpp.

78  :
79  StdExpansion(S),
80  StdExpansion1D(S),
81  StdRegions::StdSegExp(S),
82  Expansion(S),
83  Expansion1D(S),
84  m_matrixManager(S.m_matrixManager),
85  m_staticCondMatrixManager(S.m_staticCondMatrixManager)
86  {
87  }
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: SegExp.h:237
Expansion(SpatialDomains::GeometrySharedPtr pGeom)
Definition: Expansion.cpp:46
Expansion1D(SpatialDomains::Geometry1DSharedPtr pGeom)
Definition: Expansion1D.h:63
StdExpansion()
Default Constructor.
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: SegExp.h:235
Nektar::LocalRegions::SegExp::~SegExp ( )

Definition at line 93 of file SegExp.cpp.

94  {
95  }
Nektar::LocalRegions::SegExp::SegExp ( )
private

Member Function Documentation

DNekScalMatSharedPtr Nektar::LocalRegions::SegExp::CreateMatrix ( const MatrixKey mkey)
protected

Definition at line 1177 of file SegExp.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL1, ASSERTL2, Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::eFactorGaussVertex, Nektar::StdRegions::eFactorLambda, ErrorUtil::efatal, Nektar::StdRegions::eHelmholtz, Nektar::StdRegions::eHybridDGHelmBndLam, Nektar::StdRegions::eHybridDGHelmholtz, Nektar::StdRegions::eHybridDGLamToQ0, Nektar::StdRegions::eHybridDGLamToU, Nektar::StdRegions::eInterpGauss, Nektar::StdRegions::eInvHybridDGHelmholtz, Nektar::StdRegions::eInvMass, Nektar::StdRegions::eLaplacian, Nektar::StdRegions::eMass, Nektar::SpatialDomains::eNoGeomType, Nektar::StdRegions::eWeakDeriv0, Nektar::StdRegions::eWeakDeriv1, Nektar::StdRegions::eWeakDeriv2, Nektar::StdRegions::StdExpansion::GenMatrix(), Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdMatrixKey::GetConstFactors(), Nektar::StdRegions::StdMatrixKey::GetMatrixType(), Nektar::StdRegions::StdMatrixKey::GetNVarCoeff(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdMatrixKey::GetShapeType(), Nektar::StdRegions::StdExpansion::GetStdMatrix(), Nektar::StdRegions::StdMatrixKey::GetVarCoeffs(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_geom, m_matrixManager, Nektar::LocalRegions::Expansion::m_metricinfo, and NEKERROR.

1178  {
1179  DNekScalMatSharedPtr returnval;
1180  NekDouble fac;
1182 
1183  ASSERTL2(m_metricinfo->GetGtype() !=
1185  "Geometric information is not set up");
1186 
1187  switch (mkey.GetMatrixType())
1188  {
1189  case StdRegions::eMass:
1190  {
1191  if ((m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1192  || (mkey.GetNVarCoeff()))
1193  {
1194  fac = 1.0;
1195  goto UseLocRegionsMatrix;
1196  }
1197  else
1198  {
1199  fac = (m_metricinfo->GetJac(ptsKeys))[0];
1200  goto UseStdRegionsMatrix;
1201  }
1202  }
1203  break;
1204  case StdRegions::eInvMass:
1205  {
1206  if ((m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1207  || (mkey.GetNVarCoeff()))
1208  {
1209  NekDouble one = 1.0;
1210  StdRegions::StdMatrixKey masskey(
1211  StdRegions::eMass,DetShapeType(), *this);
1212  DNekMatSharedPtr mat = GenMatrix(masskey);
1213  mat->Invert();
1214 
1215  returnval = MemoryManager<DNekScalMat>::
1216  AllocateSharedPtr(one,mat);
1217  }
1218  else
1219  {
1220  fac = 1.0/(m_metricinfo->GetJac(ptsKeys))[0];
1221  goto UseStdRegionsMatrix;
1222  }
1223  }
1224  break;
1228  {
1229  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed ||
1230  mkey.GetNVarCoeff())
1231  {
1232  fac = 1.0;
1233  goto UseLocRegionsMatrix;
1234  }
1235  else
1236  {
1237  int dir = 0;
1238  switch(mkey.GetMatrixType())
1239  {
1241  dir = 0;
1242  break;
1244  ASSERTL1(m_geom->GetCoordim() >= 2,
1245  "Cannot call eWeakDeriv2 in a "
1246  "coordinate system which is not at "
1247  "least two-dimensional");
1248  dir = 1;
1249  break;
1251  ASSERTL1(m_geom->GetCoordim() == 3,
1252  "Cannot call eWeakDeriv2 in a "
1253  "coordinate system which is not "
1254  "three-dimensional");
1255  dir = 2;
1256  break;
1257  default:
1258  break;
1259  }
1260 
1261  MatrixKey deriv0key(StdRegions::eWeakDeriv0,
1262  mkey.GetShapeType(), *this);
1263 
1264  DNekMatSharedPtr WeakDerivStd = GetStdMatrix(deriv0key);
1265  fac = m_metricinfo->GetDerivFactors(ptsKeys)[dir][0]*
1266  m_metricinfo->GetJac(ptsKeys)[0];
1267 
1268  returnval = MemoryManager<DNekScalMat>::
1269  AllocateSharedPtr(fac,WeakDerivStd);
1270  }
1271  }
1272  break;
1274  {
1275  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1276  {
1277  fac = 1.0;
1278  goto UseLocRegionsMatrix;
1279  }
1280  else
1281  {
1282  int coordim = m_geom->GetCoordim();
1283  fac = 0.0;
1284  for (int i = 0; i < coordim; ++i)
1285  {
1286  fac += m_metricinfo->GetDerivFactors(ptsKeys)[i][0]*
1287  m_metricinfo->GetDerivFactors(ptsKeys)[i][0];
1288  }
1289  fac *= m_metricinfo->GetJac(ptsKeys)[0];
1290  goto UseStdRegionsMatrix;
1291  }
1292  }
1293  break;
1295  {
1296  NekDouble factor =
1297  mkey.GetConstFactor(StdRegions::eFactorLambda);
1298  MatrixKey masskey(StdRegions::eMass,
1299  mkey.GetShapeType(), *this);
1300  DNekScalMat &MassMat = *(this->m_matrixManager[masskey]);
1301  MatrixKey lapkey(StdRegions::eLaplacian, mkey.GetShapeType(),
1302  *this, mkey.GetConstFactors(),
1303  mkey.GetVarCoeffs());
1304  DNekScalMat &LapMat = *(this->m_matrixManager[lapkey]);
1305 
1306  int rows = LapMat.GetRows();
1307  int cols = LapMat.GetColumns();
1308 
1309  DNekMatSharedPtr helm =
1311 
1312  NekDouble one = 1.0;
1313  (*helm) = LapMat + factor*MassMat;
1314 
1315  returnval =
1317  }
1318  break;
1323  {
1324  NekDouble one = 1.0;
1325 
1326  DNekMatSharedPtr mat = GenMatrix(mkey);
1327  returnval =
1329  }
1330  break;
1332  {
1333  NekDouble one = 1.0;
1334 
1335 // StdRegions::StdMatrixKey hkey(StdRegions::eHybridDGHelmholtz,
1336 // DetShapeType(),*this,
1337 // mkey.GetConstant(0),
1338 // mkey.GetConstant(1));
1339  MatrixKey hkey(StdRegions::eHybridDGHelmholtz,
1340  DetShapeType(),
1341  *this, mkey.GetConstFactors(),
1342  mkey.GetVarCoeffs());
1343  DNekMatSharedPtr mat = GenMatrix(hkey);
1344 
1345  mat->Invert();
1346  returnval =
1348  }
1349  break;
1351  {
1352  DNekMatSharedPtr m_Ix;
1353  Array<OneD, NekDouble> coords(1, 0.0);
1354  StdRegions::ConstFactorMap factors = mkey.GetConstFactors();
1355  int vertex = (int)factors[StdRegions::eFactorGaussVertex];
1356 
1357  coords[0] = (vertex == 0) ? -1.0 : 1.0;
1358 
1359  m_Ix = m_base[0]->GetI(coords);
1360  returnval =
1362  }
1363  break;
1364 
1365  UseLocRegionsMatrix:
1366  {
1367  DNekMatSharedPtr mat = GenMatrix(mkey);
1368  returnval =
1370  }
1371  break;
1372  UseStdRegionsMatrix:
1373  {
1374  DNekMatSharedPtr mat = GetStdMatrix(mkey);
1375  returnval =
1377  }
1378  break;
1379  default:
1380  NEKERROR(ErrorUtil::efatal, "Matrix creation not defined");
1381  break;
1382  }
1383 
1384  return returnval;
1385  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:454
DNekMatSharedPtr GenMatrix(const StdMatrixKey &mkey)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:158
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:220
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:126
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:248
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:125
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:684
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
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
#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
Geometry is curved or has non-constant factors.
NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag > DNekScalMat
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: SegExp.h:235
DNekScalBlkMatSharedPtr Nektar::LocalRegions::SegExp::CreateStaticCondMatrix ( const MatrixKey mkey)
protected
Todo:
Really need a constructor which takes Arrays

Definition at line 1411 of file SegExp.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL2, Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::eHelmholtz, Nektar::StdRegions::eLaplacian, Nektar::SpatialDomains::eNoGeomType, Nektar::StdRegions::StdExpansion::GetBoundaryMap(), Nektar::StdRegions::StdExpansion::GetInteriorMap(), Nektar::LocalRegions::Expansion::GetLocMatrix(), Nektar::StdRegions::StdMatrixKey::GetMatrixType(), Nektar::StdRegions::StdExpansion::GetStdStaticCondMatrix(), Nektar::LocalRegions::Expansion::m_metricinfo, Nektar::StdRegions::StdExpansion::m_ncoeffs, and Nektar::StdRegions::StdExpansion::NumBndryCoeffs().

1413  {
1414  DNekScalBlkMatSharedPtr returnval;
1415 
1416  ASSERTL2(m_metricinfo->GetGtype() !=
1418  "Geometric information is not set up");
1419 
1420  // set up block matrix system
1421  int nbdry = NumBndryCoeffs();
1422  int nint = m_ncoeffs - nbdry;
1423  Array<OneD, unsigned int> exp_size(2);
1424  exp_size[0] = nbdry;
1425  exp_size[1] = nint;
1426 
1427  /// \todo Really need a constructor which takes Arrays
1429  AllocateSharedPtr(exp_size,exp_size);
1430  NekDouble factor = 1.0;
1431 
1432  switch (mkey.GetMatrixType())
1433  {
1435  case StdRegions::eHelmholtz: // special case since Helmholtz
1436  // not defined in StdRegions
1437 
1438  // use Deformed case for both regular and deformed geometries
1439  factor = 1.0;
1440  goto UseLocRegionsMatrix;
1441  break;
1442  default:
1443  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1444  {
1445  factor = 1.0;
1446  goto UseLocRegionsMatrix;
1447  }
1448  else
1449  {
1450  DNekScalMatSharedPtr mat = GetLocMatrix(mkey);
1451  factor = mat->Scale();
1452  goto UseStdRegionsMatrix;
1453  }
1454  break;
1455  UseStdRegionsMatrix:
1456  {
1457  NekDouble invfactor = 1.0/factor;
1458  NekDouble one = 1.0;
1460  DNekScalMatSharedPtr Atmp;
1461  DNekMatSharedPtr Asubmat;
1462 
1463  returnval->SetBlock(0,0,Atmp =
1464  MemoryManager<DNekScalMat>::AllocateSharedPtr(
1465  factor,Asubmat = mat->GetBlock(0,0)));
1466  returnval->SetBlock(0,1,Atmp =
1467  MemoryManager<DNekScalMat>::AllocateSharedPtr(
1468  one,Asubmat = mat->GetBlock(0,1)));
1469  returnval->SetBlock(1,0,Atmp =
1470  MemoryManager<DNekScalMat>::AllocateSharedPtr(
1471  factor,Asubmat = mat->GetBlock(1,0)));
1472  returnval->SetBlock(1,1,Atmp =
1473  MemoryManager<DNekScalMat>::AllocateSharedPtr(
1474  invfactor,Asubmat = mat->GetBlock(1,1)));
1475  }
1476  break;
1477  UseLocRegionsMatrix:
1478  {
1479  int i,j;
1480  NekDouble invfactor = 1.0/factor;
1481  NekDouble one = 1.0;
1482  DNekScalMat &mat = *GetLocMatrix(mkey);
1483  DNekMatSharedPtr A =
1485  DNekMatSharedPtr B =
1487  DNekMatSharedPtr C =
1489  DNekMatSharedPtr D =
1491 
1492  Array<OneD,unsigned int> bmap(nbdry);
1493  Array<OneD,unsigned int> imap(nint);
1494  GetBoundaryMap(bmap);
1495  GetInteriorMap(imap);
1496 
1497  for (i = 0; i < nbdry; ++i)
1498  {
1499  for (j = 0; j < nbdry; ++j)
1500  {
1501  (*A)(i,j) = mat(bmap[i],bmap[j]);
1502  }
1503 
1504  for (j = 0; j < nint; ++j)
1505  {
1506  (*B)(i,j) = mat(bmap[i],imap[j]);
1507  }
1508  }
1509 
1510  for (i = 0; i < nint; ++i)
1511  {
1512  for (j = 0; j < nbdry; ++j)
1513  {
1514  (*C)(i,j) = mat(imap[i],bmap[j]);
1515  }
1516 
1517  for (j = 0; j < nint; ++j)
1518  {
1519  (*D)(i,j) = mat(imap[i],imap[j]);
1520  }
1521  }
1522 
1523  // Calculate static condensed system
1524  if (nint)
1525  {
1526  D->Invert();
1527  (*B) = (*B)*(*D);
1528  (*A) = (*A) - (*B)*(*C);
1529  }
1530 
1531  DNekScalMatSharedPtr Atmp;
1532 
1533  returnval->SetBlock(0,0,Atmp = MemoryManager<DNekScalMat>::
1534  AllocateSharedPtr(factor,A));
1535  returnval->SetBlock(0,1,Atmp = MemoryManager<DNekScalMat>::
1536  AllocateSharedPtr(one,B));
1537  returnval->SetBlock(1,0,Atmp = MemoryManager<DNekScalMat>::
1538  AllocateSharedPtr(factor,C));
1539  returnval->SetBlock(1,1,Atmp = MemoryManager<DNekScalMat>::
1540  AllocateSharedPtr(invfactor,D));
1541  }
1542  }
1543 
1544 
1545  return returnval;
1546  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:126
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
DNekBlkMatSharedPtr GetStdStaticCondMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:689
boost::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:74
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:795
double NekDouble
boost::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:72
DNekScalMatSharedPtr GetLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:83
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:213
Geometry is curved or has non-constant factors.
NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag > DNekScalMat
void GetBoundaryMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:790
void Nektar::LocalRegions::SegExp::MultiplyByElmtInvMass ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
private
Todo:
Same method exists in ExpList and everyone references ExpList::MultiplyByElmtInvMass. Remove this one?

Definition at line 1600 of file SegExp.cpp.

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

1603  {
1604  // get Mass matrix inverse
1605  MatrixKey masskey(StdRegions::eInvMass,
1606  DetShapeType(),*this);
1607  DNekScalMatSharedPtr matsys = m_matrixManager[masskey];
1608 
1609  NekVector<NekDouble> in(m_ncoeffs,inarray,eCopy);
1610  NekVector<NekDouble> out(m_ncoeffs,outarray,eWrapper);
1611 
1612  out = (*matsys)*in;
1613  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:454
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: SegExp.h:235
void Nektar::LocalRegions::SegExp::ReverseCoeffsAndSign ( const Array< OneD, NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
private

Reverse the coefficients in a boundary interior expansion this routine is of use when we need the segment coefficients corresponding to a expansion in the reverse coordinate direction.

Definition at line 1558 of file SegExp.cpp.

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

Referenced by v_SetCoeffsToOrientation().

1561  {
1562 
1563  int m;
1564  NekDouble sgn = 1;
1565 
1566  ASSERTL1(&inarray[0] != &outarray[0],
1567  "inarray and outarray can not be the same");
1568  switch(GetBasisType(0))
1569  {
1571  //Swap vertices
1572  outarray[0] = inarray[1];
1573  outarray[1] = inarray[0];
1574  // negate odd modes
1575  for(m = 2; m < m_ncoeffs; ++m)
1576  {
1577  outarray[m] = sgn*inarray[m];
1578  sgn = -sgn;
1579  }
1580  break;
1583  for(m = 0; m < m_ncoeffs; ++m)
1584  {
1585  outarray[m_ncoeffs-1-m] = inarray[m];
1586  }
1587  break;
1588  default:
1589  ASSERTL0(false,"This basis is not allowed in this method");
1590  break;
1591  }
1592  }
#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
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
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::LocalRegions::SegExp::v_ComputeVertexNormal ( const int  vertex)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 864 of file SegExp.cpp.

References ASSERTL0, Nektar::SpatialDomains::eMovingRegular, Nektar::SpatialDomains::eRegular, Vmath::Fill(), Nektar::StdRegions::StdExpansion::GetCoordim(), Nektar::LocalRegions::Expansion::GetGeom(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion1D::m_vertexNormals, and Vmath::Smul().

865  {
866  int i;
867  const SpatialDomains::GeomFactorsSharedPtr &geomFactors =
868  GetGeom()->GetMetricInfo();
869  SpatialDomains::GeomType type = geomFactors->GetGtype();
870  const Array<TwoD, const NekDouble> &gmat =
871  geomFactors->GetDerivFactors(GetPointsKeys());
872  int nqe = m_base[0]->GetNumPoints();
873  int vCoordDim = GetCoordim();
874 
875  m_vertexNormals[vertex] =
876  Array<OneD, Array<OneD, NekDouble> >(vCoordDim);
877  Array<OneD, Array<OneD, NekDouble> > &normal =
878  m_vertexNormals[vertex];
879  for (i = 0; i < vCoordDim; ++i)
880  {
881  normal[i] = Array<OneD, NekDouble>(nqe);
882  }
883 
884  // Regular geometry case
885  if ((type == SpatialDomains::eRegular) ||
887  {
888  NekDouble vert;
889  // Set up normals
890  switch (vertex)
891  {
892  case 0:
893  for(i = 0; i < vCoordDim; ++i)
894  {
895  Vmath::Fill(nqe, -gmat[i][0], normal[i], 1);
896  }
897  break;
898  case 1:
899  for(i = 0; i < vCoordDim; ++i)
900  {
901  Vmath::Fill(nqe, gmat[i][0], normal[i], 1);
902  }
903  break;
904  default:
905  ASSERTL0(false,
906  "point is out of range (point < 2)");
907  }
908 
909  // normalise
910  vert = 0.0;
911  for (i =0 ; i < vCoordDim; ++i)
912  {
913  vert += normal[i][0]*normal[i][0];
914  }
915  vert = 1.0/sqrt(vert);
916  for (i = 0; i < vCoordDim; ++i)
917  {
918  Vmath::Smul(nqe, vert, normal[i], 1, normal[i], 1);
919  }
920  }
921  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
std::map< int, NormalVector > m_vertexNormals
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:46
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
double NekDouble
SpatialDomains::GeometrySharedPtr GetGeom() const
Definition: Expansion.cpp:148
boost::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition: GeomFactors.h:62
Geometry is straight-sided with constant geometric factors.
GeomType
Indicates the type of element geometry.
Array< OneD, LibUtilities::BasisSharedPtr > m_base
DNekMatSharedPtr Nektar::LocalRegions::SegExp::v_CreateStdMatrix ( const StdRegions::StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 1166 of file SegExp.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and Nektar::StdRegions::StdExpansion::m_base.

1168  {
1169  LibUtilities::BasisKey bkey = m_base[0]->GetBasisKey();
1172 
1173  return tmp->GetStdMatrix(mkey);
1174  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< StdSegExp > StdSegExpSharedPtr
Definition: StdSegExp.h:47
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::LocalRegions::SegExp::v_DropLocStaticCondMatrix ( const MatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1155 of file SegExp.cpp.

References m_staticCondMatrixManager.

1156  {
1157  m_staticCondMatrixManager.DeleteObject(mkey);
1158  }
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: SegExp.h:237
void Nektar::LocalRegions::SegExp::v_ExtractDataToCoeffs ( const NekDouble data,
const std::vector< unsigned int > &  nummodes,
const int  mode_offset,
NekDouble coeffs 
)
protectedvirtual

Unpack data from input file assuming it comes from.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 820 of file SegExp.cpp.

References ASSERTL0, Nektar::LibUtilities::eGauss_Lagrange, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::LibUtilities::Interp1D(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Vmath::Vcopy(), and Vmath::Zero().

825  {
826  switch(m_base[0]->GetBasisType())
827  {
829  {
830  int fillorder = min((int) nummodes[mode_offset],m_ncoeffs);
831 
832  Vmath::Zero(m_ncoeffs,coeffs,1);
833  Vmath::Vcopy(fillorder,&data[0],1,&coeffs[0],1);
834  }
835  break;
837  {
838  // Assume that input is also Gll_Lagrange
839  // but no way to check;
840  LibUtilities::PointsKey p0(
841  nummodes[mode_offset],
844  p0,data, m_base[0]->GetPointsKey(), coeffs);
845  }
846  break;
848  {
849  // Assume that input is also Gauss_Lagrange
850  // but no way to check;
851  LibUtilities::PointsKey p0(
852  nummodes[mode_offset],
855  p0,data, m_base[0]->GetPointsKey(), coeffs);
856  }
857  break;
858  default:
859  ASSERTL0(false,
860  "basis is either not set up or not hierarchicial");
861  }
862  }
#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
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:47
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
void Interp1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
this function interpolates a 1D function evaluated at the quadrature points of the basis fbasis0 to ...
Definition: Interp.cpp:54
Lagrange for SEM basis .
Definition: BasisType.h:53
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:359
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
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:50
void Nektar::LocalRegions::SegExp::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 $

Inputs:

  • inarray: array of physical quadrature points to be transformed

Outputs:

  • outarray: updated array of expansion coefficients.

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 375 of file SegExp.cpp.

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

Referenced by v_FwdTrans_BndConstrained().

378  {
379  if (m_base[0]->Collocation())
380  {
381  Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
382  }
383  else
384  {
385  v_IProductWRTBase(inarray,outarray);
386 
387  // get Mass matrix inverse
388  MatrixKey masskey(StdRegions::eInvMass, DetShapeType(), *this);
389  DNekScalMatSharedPtr matsys = m_matrixManager[masskey];
390 
391  // copy inarray in case inarray == outarray
392  NekVector<NekDouble> in(m_ncoeffs,outarray,eCopy);
393  NekVector<NekDouble> out(m_ncoeffs,outarray,eWrapper);
394 
395  out = (*matsys)*in;
396  }
397  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:454
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)->_Base[0] and return ...
Definition: SegExp.cpp:501
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
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
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: SegExp.h:235
void Nektar::LocalRegions::SegExp::v_FwdTrans_BndConstrained ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 399 of file SegExp.cpp.

References ASSERTL0, Nektar::StdRegions::StdExpansion::DetShapeType(), 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, m_staticCondMatrixManager, Nektar::StdRegions::StdExpansion::MassMatrixOp(), v_FwdTrans(), v_IProductWRTBase(), Vmath::Vcopy(), and Vmath::Vsub().

402  {
403  if(m_base[0]->Collocation())
404  {
405  Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
406  }
407  else
408  {
409  int nInteriorDofs = m_ncoeffs-2;
410  int offset;
411 
412  switch (m_base[0]->GetBasisType())
413  {
415  {
416  offset = 1;
417  }
418  break;
420  {
421  nInteriorDofs = m_ncoeffs;
422  offset = 0;
423  }
424  break;
427  {
428  offset = 2;
429  }
430  break;
431  default:
432  ASSERTL0(false,"This type of FwdTrans is not defined"
433  "for this expansion type");
434  }
435 
436  fill(outarray.get(), outarray.get()+m_ncoeffs, 0.0 );
437 
439  {
440 
441  outarray[GetVertexMap(0)] = inarray[0];
442  outarray[GetVertexMap(1)] =
443  inarray[m_base[0]->GetNumPoints()-1];
444 
445  if (m_ncoeffs>2)
446  {
447  // ideally, we would like to have tmp0 to be replaced
448  // by outarray (currently MassMatrixOp does not allow
449  // aliasing)
450  Array<OneD, NekDouble> tmp0(m_ncoeffs);
451  Array<OneD, NekDouble> tmp1(m_ncoeffs);
452 
453  StdRegions::StdMatrixKey stdmasskey(
455  MassMatrixOp(outarray,tmp0,stdmasskey);
456  v_IProductWRTBase(inarray,tmp1);
457 
458  Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
459 
460  // get Mass matrix inverse (only of interior DOF)
461  MatrixKey masskey(
463  DNekScalMatSharedPtr matsys =
464  (m_staticCondMatrixManager[masskey])->GetBlock(1,1);
465 
466  Blas::Dgemv('N',nInteriorDofs,nInteriorDofs,
467  matsys->Scale(),
468  &((matsys->GetOwnedMatrix())->GetPtr())[0],
469  nInteriorDofs,tmp1.get()+offset,1,0.0,
470  outarray.get()+offset,1);
471  }
472  }
473  else
474  {
475  SegExp::v_FwdTrans(inarray, outarray);
476 
477  }
478  }
479  }
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
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:931
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: SegExp.h:237
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: SegExp.cpp:375
Principle Modified Functions .
Definition: BasisType.h:49
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:54
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)->_Base[0] and return ...
Definition: SegExp.cpp:501
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
int GetVertexMap(const int localVertexId, bool useCoeffPacking=false)
Definition: StdExpansion.h:800
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
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
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:1038
DNekMatSharedPtr Nektar::LocalRegions::SegExp::v_GenMatrix ( const StdRegions::StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 1387 of file SegExp.cpp.

References Nektar::StdRegions::eHybridDGHelmBndLam, Nektar::StdRegions::eHybridDGHelmholtz, Nektar::StdRegions::eHybridDGLamToQ0, Nektar::StdRegions::eHybridDGLamToQ1, Nektar::StdRegions::eHybridDGLamToQ2, Nektar::StdRegions::eHybridDGLamToU, Nektar::StdRegions::StdMatrixKey::GetMatrixType(), and Nektar::LocalRegions::Expansion1D::v_GenMatrix().

1389  {
1390  DNekMatSharedPtr returnval;
1391 
1392  switch (mkey.GetMatrixType())
1393  {
1400  returnval = Expansion1D::v_GenMatrix(mkey);
1401  break;
1402  default:
1403  returnval = StdSegExp::v_GenMatrix(mkey);
1404  break;
1405  }
1406 
1407  return returnval;
1408  }
virtual DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey)
Definition: Expansion1D.cpp:43
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
const LibUtilities::BasisSharedPtr & Nektar::LocalRegions::SegExp::v_GetBasis ( int  dir) const
protectedvirtual

Definition at line 800 of file SegExp.cpp.

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

801  {
802  return GetBasis(dir);
803  }
const LibUtilities::BasisSharedPtr & GetBasis(int dir) const
This function gets the shared point to basis in the dir direction.
Definition: StdExpansion.h:118
void Nektar::LocalRegions::SegExp::v_GetCoord ( const Array< OneD, const NekDouble > &  Lcoords,
Array< OneD, NekDouble > &  coords 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 670 of file SegExp.cpp.

References ASSERTL1, and Nektar::LocalRegions::Expansion::m_geom.

673  {
674  int i;
675 
676  ASSERTL1(Lcoords[0] >= -1.0&& Lcoords[0] <= 1.0,
677  "Local coordinates are not in region [-1,1]");
678 
679  m_geom->FillGeom();
680  for(i = 0; i < m_geom->GetCoordim(); ++i)
681  {
682  coords[i] = m_geom->GetCoord(i,Lcoords);
683  }
684  }
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:125
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
int Nektar::LocalRegions::SegExp::v_GetCoordim ( void  )
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion1D.

Definition at line 774 of file SegExp.cpp.

References Nektar::LocalRegions::Expansion::m_geom.

775  {
776  return m_geom->GetCoordim();
777  }
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:125
void Nektar::LocalRegions::SegExp::v_GetCoords ( Array< OneD, NekDouble > &  coords_1,
Array< OneD, NekDouble > &  coords_2,
Array< OneD, NekDouble > &  coords_3 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 686 of file SegExp.cpp.

References Nektar::LocalRegions::Expansion::v_GetCoords().

690  {
691  Expansion::v_GetCoords(coords_0, coords_1, coords_2);
692  }
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3)
Definition: Expansion.cpp:211
DNekScalMatSharedPtr Nektar::LocalRegions::SegExp::v_GetLocMatrix ( const MatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1160 of file SegExp.cpp.

References m_matrixManager.

1161  {
1162  return m_matrixManager[mkey];
1163  }
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: SegExp.h:235
DNekScalBlkMatSharedPtr Nektar::LocalRegions::SegExp::v_GetLocStaticCondMatrix ( const MatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1149 of file SegExp.cpp.

References m_staticCondMatrixManager.

1151  {
1152  return m_staticCondMatrixManager[mkey];
1153  }
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: SegExp.h:237
int Nektar::LocalRegions::SegExp::v_GetNcoeffs ( void  ) const
protectedvirtual

Definition at line 795 of file SegExp.cpp.

References Nektar::StdRegions::StdExpansion::m_ncoeffs.

796  {
797  return m_ncoeffs;
798  }
int Nektar::LocalRegions::SegExp::v_GetNumPoints ( const int  dir) const
protectedvirtual

Definition at line 790 of file SegExp.cpp.

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

791  {
792  return GetNumPoints(dir);
793  }
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:229
const Array< OneD, const NekDouble > & Nektar::LocalRegions::SegExp::v_GetPhysNormals ( void  )
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 779 of file SegExp.cpp.

References ErrorUtil::efatal, NEKERROR, and Nektar::NullNekDouble1DArray.

780  {
781  NEKERROR(ErrorUtil::efatal, "Got to SegExp");
782  return NullNekDouble1DArray;
783  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:158
static Array< OneD, NekDouble > NullNekDouble1DArray
StdRegions::Orientation Nektar::LocalRegions::SegExp::v_GetPorient ( int  point)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 762 of file SegExp.cpp.

References Nektar::LocalRegions::Expansion::m_geom.

763  {
764  return m_geom->GetPorient(point);
765  }
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:125
StdRegions::StdExpansionSharedPtr Nektar::LocalRegions::SegExp::v_GetStdExp ( void  ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 768 of file SegExp.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and Nektar::StdRegions::StdExpansion::m_base.

769  {
771  ::AllocateSharedPtr(m_base[0]->GetBasisKey());
772  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::LocalRegions::SegExp::v_GetVertexPhysVals ( const int  vertex,
const Array< OneD, const NekDouble > &  inarray,
NekDouble outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 695 of file SegExp.cpp.

References Vmath::Ddot(), Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::StdRegions::eFactorGaussVertex, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::StdRegions::eInterpGauss, Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::StdRegions::StdExpansion::m_base, and m_matrixManager.

699  {
700  int nquad = m_base[0]->GetNumPoints();
701 
703  {
704  switch (vertex)
705  {
706  case 0:
707  outarray = inarray[0];
708  break;
709  case 1:
710  outarray = inarray[nquad - 1];
711  break;
712  }
713  }
714  else
715  {
717  factors[StdRegions::eFactorGaussVertex] = vertex;
718 
719  StdRegions::StdMatrixKey key(
721  DetShapeType(),*this,factors);
722 
723  DNekScalMatSharedPtr mat_gauss = m_matrixManager[key];
724 
725  outarray = Blas::Ddot(nquad, mat_gauss->GetOwnedMatrix()
726  ->GetPtr().get(), 1, &inarray[0], 1);
727  }
728  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:454
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:248
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:47
T Ddot(int n, const Array< OneD, const T > &w, const int incw, const Array< OneD, const T > &x, const int incx, const Array< OneD, const int > &y, const int incy)
Definition: VmathArray.hpp:434
LibUtilities::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
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: SegExp.h:235
void Nektar::LocalRegions::SegExp::v_HelmholtzMatrixOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdRegions::StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 1034 of file SegExp.cpp.

References ASSERTL0, Nektar::StdRegions::StdExpansion::BwdTrans(), Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_geom, Nektar::LocalRegions::Expansion::m_metricinfo, Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::StdRegions::StdExpansion::PhysDeriv(), v_IProductWRTBase(), Vmath::Vmul(), and Vmath::Vvtvp().

1038  {
1039  int nquad = m_base[0]->GetNumPoints();
1040  const Array<TwoD, const NekDouble>& gmat =
1041  m_metricinfo->GetDerivFactors(GetPointsKeys());
1042  const NekDouble lambda =
1043  mkey.GetConstFactor(StdRegions::eFactorLambda);
1044 
1045  Array<OneD,NekDouble> physValues(nquad);
1046  Array<OneD,NekDouble> dPhysValuesdx(nquad);
1047  Array<OneD,NekDouble> wsp(m_ncoeffs);
1048 
1049  BwdTrans(inarray, physValues);
1050 
1051  // mass matrix operation
1052  v_IProductWRTBase((m_base[0]->GetBdata()),physValues,wsp,1);
1053 
1054  // Laplacian matrix operation
1055  switch (m_geom->GetCoordim())
1056  {
1057  case 1:
1058  {
1059  PhysDeriv(physValues,dPhysValuesdx);
1060 
1061  // multiply with the proper geometric factors
1062  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1063  {
1064  Vmath::Vmul(nquad,
1065  &gmat[0][0],1,dPhysValuesdx.get(),1,
1066  dPhysValuesdx.get(),1);
1067  }
1068  else
1069  {
1070  Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.get(), 1);
1071  }
1072  }
1073  break;
1074  case 2:
1075  {
1076  Array<OneD,NekDouble> dPhysValuesdy(nquad);
1077 
1078  PhysDeriv(physValues, dPhysValuesdx, dPhysValuesdy);
1079 
1080  // multiply with the proper geometric factors
1081  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1082  {
1083  Vmath::Vmul (nquad,
1084  &gmat[0][0], 1, dPhysValuesdx.get(), 1,
1085  dPhysValuesdx.get(), 1);
1086  Vmath::Vvtvp(nquad,
1087  &gmat[1][0], 1, dPhysValuesdy.get(), 1,
1088  dPhysValuesdx.get(), 1,
1089  dPhysValuesdx.get(), 1);
1090  }
1091  else
1092  {
1093  Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.get(), 1);
1094  Blas::Daxpy(nquad,
1095  gmat[1][0], dPhysValuesdy.get(), 1,
1096  dPhysValuesdx.get(), 1);
1097  }
1098  }
1099  break;
1100  case 3:
1101  {
1102  Array<OneD,NekDouble> dPhysValuesdy(nquad);
1103  Array<OneD,NekDouble> dPhysValuesdz(nquad);
1104 
1105  PhysDeriv(physValues, dPhysValuesdx,
1106  dPhysValuesdy, dPhysValuesdz);
1107 
1108  // multiply with the proper geometric factors
1109  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1110  {
1111  Vmath::Vmul (nquad,
1112  &gmat[0][0], 1, dPhysValuesdx.get(), 1,
1113  dPhysValuesdx.get(), 1);
1114  Vmath::Vvtvp(nquad,
1115  &gmat[1][0], 1, dPhysValuesdy.get(), 1,
1116  dPhysValuesdx.get(), 1,
1117  dPhysValuesdx.get(), 1);
1118  Vmath::Vvtvp(nquad,
1119  &gmat[2][0], 1, dPhysValuesdz.get(), 1,
1120  dPhysValuesdx.get(), 1,
1121  dPhysValuesdx.get(), 1);
1122  }
1123  else
1124  {
1125  Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.get(), 1);
1126  Blas::Daxpy(nquad,
1127  gmat[1][0], dPhysValuesdy.get(), 1,
1128  dPhysValuesdx.get(), 1);
1129  Blas::Daxpy(nquad,
1130  gmat[2][0], dPhysValuesdz.get(),
1131  1, dPhysValuesdx.get(), 1);
1132  }
1133  }
1134  break;
1135  default:
1136  ASSERTL0(false,"Wrong number of dimensions");
1137  break;
1138  }
1139 
1140  v_IProductWRTBase(m_base[0]->GetDbdata(),dPhysValuesdx,outarray,1);
1141  Blas::Daxpy(m_ncoeffs, lambda, wsp.get(), 1, outarray.get(), 1);
1142  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:428
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:126
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:125
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)->_Base[0] and return ...
Definition: SegExp.cpp:501
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)
double NekDouble
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Backward transformation from coefficient space to physical space...
Definition: StdExpansion.h:509
Array< OneD, LibUtilities::BasisSharedPtr > m_base
Geometry is curved or has non-constant factors.
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
NekDouble Nektar::LocalRegions::SegExp::v_Integral ( const Array< OneD, const NekDouble > &  inarray)
protectedvirtual

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

Inputs:

  • inarray: definition of function to be returned at quadrature point of expansion.

Outputs:

  • returns $\int^1_{-1} u(\xi_1)d \xi_1 $ where $inarray[i] = u(\xi_{1i}) $

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 116 of file SegExp.cpp.

References Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_metricinfo, Vmath::Smul(), and Vmath::Vmul().

118  {
119  int nquad0 = m_base[0]->GetNumPoints();
120  Array<OneD, const NekDouble> jac = m_metricinfo->GetJac(GetPointsKeys());
121  NekDouble ival;
122  Array<OneD,NekDouble> tmp(nquad0);
123 
124  // multiply inarray with Jacobian
125  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
126  {
127  Vmath::Vmul(nquad0, jac, 1, inarray, 1, tmp,1);
128  }
129  else
130  {
131  Vmath::Smul(nquad0, jac[0], inarray, 1, tmp, 1);
132  }
133 
134  // call StdSegExp version;
135  ival = StdSegExp::v_Integral(tmp);
136  //ival = StdSegExp::Integral(tmp);
137  return ival;
138  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:126
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
double NekDouble
Array< OneD, LibUtilities::BasisSharedPtr > m_base
Geometry is curved or has non-constant factors.
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::LocalRegions::SegExp::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)->_Base[0] and return in outarray.

Wrapper call to SegExp::IProduct_WRT_B

Input:

  • inarray: array of function evaluated at the physical collocation points

Output:

  • outarray: array of inner product with respect to each basis over region

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 501 of file SegExp.cpp.

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

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

504  {
505  v_IProductWRTBase(m_base[0]->GetBdata(),inarray,outarray,1);
506  }
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)->_Base[0] and return ...
Definition: SegExp.cpp:501
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Nektar::LocalRegions::SegExp::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}) $.

Inputs:

  • base: an array definiing the local basis for the inner product usually passed from Basis->get_bdata() or Basis->get_Dbdata()
  • inarray: physical point array of function to be integrated $ u(\xi_1) $
  • coll_check: Flag 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)

Output:

  • outarray: array of coefficients representing the inner product of function with ever mode in the exapnsion

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 536 of file SegExp.cpp.

References Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_metricinfo, Vmath::Smul(), and Vmath::Vmul().

541  {
542  int nquad0 = m_base[0]->GetNumPoints();
543  Array<OneD, const NekDouble> jac = m_metricinfo->GetJac(GetPointsKeys());
544  Array<OneD,NekDouble> tmp(nquad0);
545 
546 
547  // multiply inarray with Jacobian
548  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
549  {
550  Vmath::Vmul(nquad0, jac, 1, inarray, 1, tmp, 1);
551  }
552  else
553  {
554  Vmath::Smul(nquad0, jac[0], inarray, 1, tmp, 1);
555  }
556  StdSegExp::v_IProductWRTBase(base,tmp,outarray,coll_check);
557  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:126
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
Array< OneD, LibUtilities::BasisSharedPtr > m_base
Geometry is curved or has non-constant factors.
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::LocalRegions::SegExp::v_IProductWRTDerivBase ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 560 of file SegExp.cpp.

References ASSERTL1, Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_geom, Nektar::LocalRegions::Expansion::m_metricinfo, Vmath::Smul(), v_IProductWRTBase(), and Vmath::Vmul().

564  {
565  int nquad = m_base[0]->GetNumPoints();
566  const Array<TwoD, const NekDouble>& gmat =
567  m_metricinfo->GetDerivFactors(GetPointsKeys());
568 
569  Array<OneD, NekDouble> tmp1(nquad);
570 
571  switch(dir)
572  {
573  case 0:
574  {
575  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
576  {
577  Vmath::Vmul(nquad,gmat[0],1,inarray,1,tmp1,1);
578  }
579  else
580  {
581  Vmath::Smul(nquad, gmat[0][0], inarray, 1, tmp1, 1);
582  }
583  }
584  break;
585  case 1:
586  {
587  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
588  {
589  Vmath::Vmul(nquad,gmat[1],1,inarray,1,tmp1,1);
590  }
591  else
592  {
593  Vmath::Smul(nquad, gmat[1][0], inarray, 1, tmp1, 1);
594  }
595  }
596  break;
597  case 2:
598  {
599  ASSERTL1(m_geom->GetCoordim() == 3,"input dir is out of range");
600  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
601  {
602  Vmath::Vmul(nquad,gmat[2],1,inarray,1,tmp1,1);
603  }
604  else
605  {
606  Vmath::Smul(nquad, gmat[2][0], inarray, 1, tmp1, 1);
607  }
608  }
609  break;
610  default:
611  {
612  ASSERTL1(false,"input dir is out of range");
613  }
614  break;
615  }
616  v_IProductWRTBase(m_base[0]->GetDbdata(),tmp1,outarray,1);
617  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:126
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:125
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)->_Base[0] and return ...
Definition: SegExp.cpp:501
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
#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
Geometry is curved or has non-constant factors.
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::LocalRegions::SegExp::v_LaplacianMatrixOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdRegions::StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 928 of file SegExp.cpp.

References ASSERTL0, Nektar::StdRegions::StdExpansion::BwdTrans(), Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_geom, Nektar::LocalRegions::Expansion::m_metricinfo, Nektar::StdRegions::StdExpansion::PhysDeriv(), v_IProductWRTBase(), Vmath::Vmul(), and Vmath::Vvtvp().

932  {
933  int nquad = m_base[0]->GetNumPoints();
934  const Array<TwoD, const NekDouble>& gmat =
935  m_metricinfo->GetDerivFactors(GetPointsKeys());
936 
937  Array<OneD,NekDouble> physValues(nquad);
938  Array<OneD,NekDouble> dPhysValuesdx(nquad);
939 
940  BwdTrans(inarray,physValues);
941 
942  // Laplacian matrix operation
943  switch (m_geom->GetCoordim())
944  {
945  case 1:
946  {
947  PhysDeriv(physValues,dPhysValuesdx);
948 
949  // multiply with the proper geometric factors
950  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
951  {
952  Vmath::Vmul(nquad,
953  &gmat[0][0],1,dPhysValuesdx.get(),1,
954  dPhysValuesdx.get(),1);
955  }
956  else
957  {
958  Blas::Dscal(nquad,
959  gmat[0][0], dPhysValuesdx.get(), 1);
960  }
961  }
962  break;
963  case 2:
964  {
965  Array<OneD,NekDouble> dPhysValuesdy(nquad);
966 
967  PhysDeriv(physValues,dPhysValuesdx,dPhysValuesdy);
968 
969  // multiply with the proper geometric factors
970  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
971  {
972  Vmath::Vmul (nquad,
973  &gmat[0][0],1,dPhysValuesdx.get(),1,
974  dPhysValuesdx.get(),1);
975  Vmath::Vvtvp(nquad,
976  &gmat[1][0],1,dPhysValuesdy.get(),1,
977  dPhysValuesdx.get(),1,
978  dPhysValuesdx.get(),1);
979  }
980  else
981  {
982  Blas::Dscal(nquad,
983  gmat[0][0], dPhysValuesdx.get(), 1);
984  Blas::Daxpy(nquad,
985  gmat[1][0], dPhysValuesdy.get(), 1,
986  dPhysValuesdx.get(), 1);
987  }
988  }
989  break;
990  case 3:
991  {
992  Array<OneD,NekDouble> dPhysValuesdy(nquad);
993  Array<OneD,NekDouble> dPhysValuesdz(nquad);
994 
995  PhysDeriv(physValues,dPhysValuesdx,
996  dPhysValuesdy,dPhysValuesdz);
997 
998  // multiply with the proper geometric factors
999  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1000  {
1001  Vmath::Vmul (nquad,
1002  &gmat[0][0], 1, dPhysValuesdx.get(), 1,
1003  dPhysValuesdx.get(),1);
1004  Vmath::Vvtvp(nquad,
1005  &gmat[1][0], 1, dPhysValuesdy.get(), 1,
1006  dPhysValuesdx.get(),1,
1007  dPhysValuesdx.get(),1);
1008  Vmath::Vvtvp(nquad,
1009  &gmat[2][0], 1, dPhysValuesdz.get(), 1,
1010  dPhysValuesdx.get(),1,
1011  dPhysValuesdx.get(),1);
1012  }
1013  else
1014  {
1015  Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.get(), 1);
1016  Blas::Daxpy(nquad,
1017  gmat[1][0], dPhysValuesdy.get(), 1,
1018  dPhysValuesdx.get(), 1);
1019  Blas::Daxpy(nquad,
1020  gmat[2][0], dPhysValuesdz.get(), 1,
1021  dPhysValuesdx.get(), 1);
1022  }
1023  }
1024  break;
1025  default:
1026  ASSERTL0(false,"Wrong number of dimensions");
1027  break;
1028  }
1029 
1030  v_IProductWRTBase(m_base[0]->GetDbdata(),dPhysValuesdx,outarray,1);
1031  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:428
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:126
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:125
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)->_Base[0] and return ...
Definition: SegExp.cpp:501
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 BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Backward transformation from coefficient space to physical space...
Definition: StdExpansion.h:509
Array< OneD, LibUtilities::BasisSharedPtr > m_base
Geometry is curved or has non-constant factors.
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
SpatialDomains::GeomType Nektar::LocalRegions::SegExp::v_MetricInfoType ( )
protectedvirtual

Definition at line 785 of file SegExp.cpp.

References Nektar::LocalRegions::Expansion::m_metricinfo.

786  {
787  return m_metricinfo->GetGtype();
788  }
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:126
void Nektar::LocalRegions::SegExp::v_NormVectorIProductWRTBase ( const Array< OneD, const NekDouble > &  Fx,
const Array< OneD, const NekDouble > &  Fy,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 619 of file SegExp.cpp.

References Nektar::StdRegions::StdExpansion::GetEdgeNormal(), Nektar::LocalRegions::Expansion1D::GetLeftAdjacentElementEdge(), Nektar::LocalRegions::Expansion1D::GetLeftAdjacentElementExp(), Nektar::StdRegions::StdExpansion::m_base, v_IProductWRTBase(), Vmath::Vmul(), and Vmath::Vvtvp().

623  {
624  int nq = m_base[0]->GetNumPoints();
625  Array<OneD, NekDouble > Fn(nq);
626 // cout << "I am segment " << GetGeom()->GetGlobalID() << endl;
627 // cout << "I want edge " << GetLeftAdjacentElementEdge() << endl;
628 // @TODO: This routine no longer makes sense as a normal is not unique to an edge
629  const Array<OneD, const Array<OneD, NekDouble> >
630  &normals =
633  Vmath::Vmul (nq, &Fx[0], 1, &normals[0][0], 1, &Fn[0], 1);
634  Vmath::Vvtvp(nq, &Fy[0], 1, &normals[1][0], 1, &Fn[0], 1, &Fn[0], 1);
635 
636  v_IProductWRTBase(Fn,outarray);
637  }
const NormalVector & GetEdgeNormal(const int edge) const
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:428
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)->_Base[0] and return ...
Definition: SegExp.cpp:501
Expansion2DSharedPtr GetLeftAdjacentElementExp() const
Definition: Expansion1D.h:116
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::LocalRegions::SegExp::v_NumBndryCoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 806 of file SegExp.cpp.

807  {
808  return 2;
809  }
int Nektar::LocalRegions::SegExp::v_NumDGBndryCoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 812 of file SegExp.cpp.

813  {
814  return 2;
815  }
void Nektar::LocalRegions::SegExp::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

Input:

  • n: number of derivatives to be evaluated where $ n \leq dim$
  • inarray: array of function evaluated at the quadrature points

Output:

  • outarray: array of the derivatives $ du/d_{\xi_1}|_{\xi_{1i}} d\xi_1/dx, du/d_{\xi_1}|_{\xi_{1i}} d\xi_1/dy, du/d_{\xi_1}|_{\xi_{1i}} d\xi_1/dz, $ depending on value of dim

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 164 of file SegExp.cpp.

References Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_metricinfo, Nektar::StdRegions::StdExpansion1D::PhysTensorDeriv(), Vmath::Smul(), and Vmath::Vmul().

Referenced by v_PhysDeriv_n().

169  {
170  int nquad0 = m_base[0]->GetNumPoints();
171  Array<TwoD, const NekDouble> gmat =
172  m_metricinfo->GetDerivFactors(GetPointsKeys());
173  Array<OneD,NekDouble> diff(nquad0);
174 
175  //StdExpansion1D::PhysTensorDeriv(inarray,diff);
176  PhysTensorDeriv(inarray,diff);
177  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
178  {
179  if(out_d0.num_elements())
180  {
181  Vmath::Vmul(nquad0,&gmat[0][0],1,&diff[0],1,
182  &out_d0[0],1);
183  }
184 
185  if(out_d1.num_elements())
186  {
187  Vmath::Vmul(nquad0,&gmat[1][0],1,&diff[0],1,
188  &out_d1[0],1);
189  }
190 
191  if(out_d2.num_elements())
192  {
193  Vmath::Vmul(nquad0,&gmat[2][0],1,&diff[0],1,
194  &out_d2[0],1);
195  }
196  }
197  else
198  {
199  if(out_d0.num_elements())
200  {
201  Vmath::Smul(nquad0, gmat[0][0], diff, 1,
202  out_d0, 1);
203  }
204 
205  if(out_d1.num_elements())
206  {
207  Vmath::Smul(nquad0, gmat[1][0], diff, 1,
208  out_d1, 1);
209  }
210 
211  if(out_d2.num_elements())
212  {
213  Vmath::Smul(nquad0, gmat[2][0], diff, 1,
214  out_d2, 1);
215  }
216  }
217  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:126
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
void 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...
Array< OneD, LibUtilities::BasisSharedPtr > m_base
Geometry is curved or has non-constant factors.
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::LocalRegions::SegExp::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::StdSegExp.

Definition at line 317 of file SegExp.cpp.

References ASSERTL1, Nektar::NullNekDouble1DArray, and Nektar::StdRegions::StdExpansion::PhysDeriv().

320  {
321  switch(dir)
322  {
323  case 0:
324  {
325  PhysDeriv(inarray, outarray, NullNekDouble1DArray,
327  }
328  break;
329  case 1:
330  {
331  PhysDeriv(inarray, NullNekDouble1DArray, outarray,
333  }
334  break;
335  case 2:
336  {
338  NullNekDouble1DArray, outarray);
339  }
340  break;
341  default:
342  {
343  ASSERTL1(false,"input dir is out of range");
344  }
345  break;
346  }
347  }
static Array< OneD, NekDouble > NullNekDouble1DArray
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)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
void Nektar::LocalRegions::SegExp::v_PhysDeriv_n ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_dn 
)
protectedvirtual

Evaluate the derivative normal to a line: $ d/dn=\frac{spacedim}{||normal||}d/d{\xi} $. The derivative is calculated performing the product $ du/d{s}=\nabla u \cdot normal $.

Parameters
inarrayfunction to derive
out_dnresult of the derivative operation

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 267 of file SegExp.cpp.

References ASSERTL0, Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_geom, Nektar::LocalRegions::Expansion::m_metricinfo, Nektar::NullNekDoubleArrayofArray, v_PhysDeriv(), Vmath::Vadd(), Vmath::Vmul(), and Vmath::Zero().

270  {
271  int nquad0 = m_base[0]->GetNumPoints();
272  Array<TwoD, const NekDouble> gmat =
273  m_metricinfo->GetDerivFactors(GetPointsKeys());
274  int coordim = m_geom->GetCoordim();
275  Array<OneD, NekDouble> out_dn_tmp(nquad0,0.0);
276  switch(coordim)
277  {
278  case 2:
279 
280  Array<OneD, NekDouble> inarray_d0(nquad0);
281  Array<OneD, NekDouble> inarray_d1(nquad0);
282 
283  v_PhysDeriv(inarray,inarray_d0,inarray_d1);
284  Array<OneD, Array<OneD, NekDouble> > normals;
285  normals = Array<OneD, Array<OneD, NekDouble> >(coordim);
286  cout<<"der_n"<<endl;
287  for(int k=0; k<coordim; ++k)
288  {
289  normals[k]= Array<OneD, NekDouble>(nquad0);
290  }
291 // @TODO: this routine no longer makes sense, since normals are not unique on
292 // an edge
293 // normals = GetMetricInfo()->GetNormal();
294  for(int i=0; i<nquad0; i++)
295  {
296 cout<<"nx= "<<normals[0][i]<<" ny="<<normals[1][i]<<endl;
297  }
299  "normal vectors do not exist: check if a"
300  "boundary region is defined as I ");
301  // \nabla u \cdot normal
302  Vmath::Vmul(nquad0,normals[0],1,inarray_d0,1,out_dn_tmp,1);
303  Vmath::Vadd(nquad0,out_dn_tmp,1,out_dn,1,out_dn,1);
304  Vmath::Zero(nquad0,out_dn_tmp,1);
305  Vmath::Vmul(nquad0,normals[1],1,inarray_d1,1,out_dn_tmp,1);
306  Vmath::Vadd(nquad0,out_dn_tmp,1,out_dn,1,out_dn,1);
307 
308  for(int i=0; i<nquad0; i++)
309  {
310 cout<<"deps/dx ="<<inarray_d0[i]<<" deps/dy="<<inarray_d1[i]<<endl;
311  }
312 
313 
314  }
315 
316  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:126
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:125
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: SegExp.cpp:164
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayofArray
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:359
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:285
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:169
void Nektar::LocalRegions::SegExp::v_PhysDeriv_s ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_ds 
)
protectedvirtual

Evaluate the derivative along a line: $ d/ds=\frac{spacedim}{||tangent||}d/d{\xi} $. The derivative is calculated performing the product $ du/d{s}=\nabla u \cdot tangent $.

Parameters
inarrayfunction to derive
out_dsresult of the derivative operation

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 227 of file SegExp.cpp.

References Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_geom, Nektar::LocalRegions::Expansion::m_metricinfo, Nektar::StdRegions::StdExpansion1D::PhysTensorDeriv(), Vmath::Smul(), Vmath::Vdiv(), and Vmath::Zero().

230  {
231  int nquad0 = m_base[0]->GetNumPoints();
232  int coordim = m_geom->GetCoordim();
233  Array<OneD, NekDouble> diff (nquad0);
234  //this operation is needed if you put out_ds==inarray
235  Vmath::Zero(nquad0,out_ds,1);
236  switch(coordim)
237  {
238  case 2:
239  //diff= dU/de
240  Array<OneD,NekDouble> diff(nquad0);
241 
242  PhysTensorDeriv(inarray,diff);
243 
244  //get dS/de= (Jac)^-1
245  Array<OneD, NekDouble> Jac = m_metricinfo->GetJac(GetPointsKeys());
246  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
247  {
248  //calculate the derivative as (dU/de)*(Jac)^-1
249  Vmath::Vdiv(nquad0,diff,1,Jac ,1,out_ds,1);
250  }
251  else
252  {
253  NekDouble invJac = 1/Jac[0];
254  Vmath::Smul(nquad0, invJac,diff,1,out_ds,1);
255  }
256  }
257  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:126
void Vdiv(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:227
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:125
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
double NekDouble
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 Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:359
Array< OneD, LibUtilities::BasisSharedPtr > m_base
Geometry is curved or has non-constant factors.
NekDouble Nektar::LocalRegions::SegExp::v_PhysEvaluate ( const Array< OneD, const NekDouble > &  coord,
const Array< OneD, const NekDouble > &  physvals 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 657 of file SegExp.cpp.

References ASSERTL0, and Nektar::LocalRegions::Expansion::m_geom.

660  {
661  Array<OneD,NekDouble> Lcoord = Array<OneD,NekDouble>(1);
662 
663  ASSERTL0(m_geom,"m_geom not defined");
664  m_geom->GetLocCoords(coord,Lcoord);
665 
666  return StdSegExp::v_PhysEvaluate(Lcoord, physvals);
667  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:125
void Nektar::LocalRegions::SegExp::v_SetCoeffsToOrientation ( Array< OneD, NekDouble > &  coeffs,
StdRegions::Orientation  dir 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 735 of file SegExp.cpp.

738  {
739  v_SetCoeffsToOrientation(dir,coeffs,coeffs);
740  }
virtual void v_SetCoeffsToOrientation(Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
Definition: SegExp.cpp:735
void Nektar::LocalRegions::SegExp::v_SetCoeffsToOrientation ( StdRegions::Orientation  dir,
Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 742 of file SegExp.cpp.

References Nektar::StdRegions::eBackwards, and ReverseCoeffsAndSign().

746  {
747 
748  if (dir == StdRegions::eBackwards)
749  {
750  if (&inarray[0] != &outarray[0])
751  {
752  Array<OneD,NekDouble> intmp (inarray);
753  ReverseCoeffsAndSign(intmp,outarray);
754  }
755  else
756  {
757  ReverseCoeffsAndSign(inarray,outarray);
758  }
759  }
760  }
void ReverseCoeffsAndSign(const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Reverse the coefficients in a boundary interior expansion this routine is of use when we need the seg...
Definition: SegExp.cpp:1558
NekDouble Nektar::LocalRegions::SegExp::v_StdPhysEvaluate ( const Array< OneD, const NekDouble > &  Lcoord,
const Array< OneD, const NekDouble > &  physvals 
)
protectedvirtual

Given the local cartesian coordinate Lcoord evaluate the value of physvals at this point by calling through to the StdExpansion method

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 649 of file SegExp.cpp.

652  {
653  // Evaluate point in local (eta) coordinates.
654  return StdSegExp::v_PhysEvaluate(Lcoord,physvals);
655  }

Member Data Documentation

LibUtilities::NekManager<MatrixKey, DNekScalMat, MatrixKey::opLess> Nektar::LocalRegions::SegExp::m_matrixManager
private
LibUtilities::NekManager<MatrixKey, DNekScalBlkMat, MatrixKey::opLess> Nektar::LocalRegions::SegExp::m_staticCondMatrixManager
private