Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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::BasisSharedPtr
GetBasis (int dir) const
 This function gets the shared point to basis in the dir direction. More...
 
int GetNcoeffs (void) const
 This function returns the total number of coefficients used in the expansion. More...
 
int GetTotPoints () const
 This function returns the total number of quadrature points used in the element. More...
 
LibUtilities::BasisType GetBasisType (const int dir) const
 This function returns the type of basis used in the dir direction. More...
 
int GetBasisNumModes (const int dir) const
 This function returns the number of expansion modes in the dir direction. More...
 
int EvalBasisNumModesMax (void) const
 This function returns the maximum number of expansion modes over all local directions. More...
 
LibUtilities::PointsType GetPointsType (const int dir) const
 This function returns the type of quadrature points used in the dir direction. More...
 
int GetNumPoints (const int dir) const
 This function returns the number of quadrature points in the dir direction. More...
 
const Array< OneD, const
NekDouble > & 
GetPoints (const int dir) const
 This function returns a pointer to the array containing the quadrature points in dir direction. More...
 
int GetNverts () const
 This function returns the number of vertices of the expansion domain. More...
 
int GetNedges () const
 This function returns the number of edges of the expansion domain. More...
 
int GetEdgeNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th edge. More...
 
int GetTotalEdgeIntNcoeffs () const
 
int GetEdgeNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th edge. More...
 
int DetCartesianDirOfEdge (const int edge)
 
const LibUtilities::BasisKey DetEdgeBasisKey (const int i) const
 
const LibUtilities::BasisKey DetFaceBasisKey (const int i, const int k) const
 
int GetFaceNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th face. More...
 
int GetFaceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th face. More...
 
int GetFaceIntNcoeffs (const int i) const
 
int GetTotalFaceIntNcoeffs () const
 
int GetTraceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th edge/face. More...
 
LibUtilities::PointsKey GetFacePointsKey (const int i, const int j) const
 
int NumBndryCoeffs (void) const
 
int NumDGBndryCoeffs (void) const
 
LibUtilities::BasisType GetEdgeBasisType (const int i) const
 This function returns the type of expansion basis on the i-th edge. More...
 
const LibUtilities::PointsKey GetNodalPointsKey () const
 This function returns the type of expansion Nodal point type if defined. More...
 
int GetNfaces () const
 This function returns the number of faces of the expansion domain. More...
 
int GetNtrace () const
 Returns the number of trace elements connected to this element. More...
 
LibUtilities::ShapeType DetShapeType () const
 This function returns the shape of the expansion domain. More...
 
boost::shared_ptr< StdExpansionGetStdExp (void) const
 
boost::shared_ptr< StdExpansionGetLinStdExp (void) const
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion ()
 
bool IsNodalNonTensorialExp ()
 
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space. More...
 
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space. More...
 
void FwdTrans_BndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble Integral (const Array< OneD, const NekDouble > &inarray)
 This function integrates the specified function over the domain. More...
 
void FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 This function fills the array outarray with the mode-th mode of the expansion. More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 this function calculates the inner product of a given function f with the different modes of the expansion More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check)
 
void IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
int GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void SetElmtId (const int id)
 Set the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2=NullNekDouble1DArray, Array< OneD, NekDouble > &coords_3=NullNekDouble1DArray)
 this function returns the physical coordinates of the quadrature points of the expansion More...
 
void GetCoord (const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
 given the coordinates of a point of the element in the local collapsed coordinate system, this function calculates the physical coordinates of the point More...
 
DNekMatSharedPtr GetStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr GetStdStaticCondMatrix (const StdMatrixKey &mkey)
 
IndexMapValuesSharedPtr GetIndexMap (const IndexMapKey &ikey)
 
const Array< OneD, const
NekDouble > & 
GetPhysNormals (void)
 
void SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
virtual void SetUpPhysNormals (const int edge)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
void DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
StdRegions::Orientation GetForient (int face)
 
StdRegions::Orientation GetEorient (int edge)
 
StdRegions::Orientation GetPorient (int point)
 
StdRegions::Orientation GetCartesianEorient (int edge)
 
void SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
 
void SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
int CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
NekDouble StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
int GetCoordim ()
 
void GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
void GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
int GetVertexMap (const int localVertexId, bool useCoeffPacking=false)
 
void GetEdgeInteriorMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
 
void GetFaceNumModes (const int fid, const Orientation faceOrient, int &numModes0, int &numModes1)
 
void GetFaceInteriorMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
 
void GetEdgeToElementMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int P=-1)
 
void GetFaceToElementMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int nummodesA=-1, int nummodesB=-1)
 
void GetEdgePhysVals (const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Extract the physical values along edge edge from inarray into outarray following the local edge orientation and point distribution defined by defined in EdgeExp. More...
 
void GetEdgePhysVals (const int edge, const boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetTracePhysVals (const int edge, const boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetVertexPhysVals (const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray)
 
void GetEdgeInterpVals (const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetEdgeQFactors (const int edge, Array< OneD, NekDouble > &outarray)
 Extract the metric factors to compute the contravariant fluxes along edge edge and stores them into outarray following the local edge orientation (i.e. anticlockwise convention). More...
 
void GetFacePhysVals (const int face, const boost::shared_ptr< StdExpansion > &FaceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=eNoOrientation)
 
void GetEdgePhysMap (const int edge, Array< OneD, int > &outarray)
 
void GetFacePhysMap (const int face, Array< OneD, int > &outarray)
 
void MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
DNekMatSharedPtr CreateGeneralMatrix (const StdMatrixKey &mkey)
 this function generates the mass matrix $\mathbf{M}[i][j] = \int \phi_i(\mathbf{x}) \phi_j(\mathbf{x}) d\mathbf{x}$ More...
 
void GeneralMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
DNekMatSharedPtr GenMatrix (const StdMatrixKey &mkey)
 
void PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void PhysDeriv_s (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds)
 
void PhysDeriv_n (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn)
 
void PhysDirectionalDeriv (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &outarray)
 
void StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void AddRobinMassMatrix (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
void AddRobinEdgeContribution (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, Array< OneD, NekDouble > &coeffs)
 
NekDouble PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta. More...
 
const boost::shared_ptr
< SpatialDomains::GeomFactors > & 
GetMetricInfo (void) const
 
virtual int v_GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
virtual void v_SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, 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)
 
bool FaceNormalNegated (const int face)
 
void ComputeVertexNormal (const int vertex)
 
const NormalVectorGetFaceNormal (const int face) const
 
const NormalVectorGetVertexNormal (const int vertex) const
 
const NormalVectorGetSurfaceNormal (const int id) const
 
const LibUtilities::PointsKeyVector GetPointsKeys () const
 
Array< OneD, unsigned int > GetEdgeInverseBoundaryMap (int eid)
 
Array< OneD, unsigned int > GetFaceInverseBoundaryMap (int fid, StdRegions::Orientation faceOrient=eNoOrientation)
 
DNekMatSharedPtr BuildInverseTransformationMatrix (const DNekScalMatSharedPtr &m_transformationmatrix)
 
void PhysInterpToSimplexEquiSpaced (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int npset=-1)
 This function performs an interpolation from the physical space points provided at input into an array of equispaced points which are not the collapsed coordinate. So for a tetrahedron you will only get a tetrahedral number of values. More...
 
void GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true)
 This function provides the connectivity of local simplices (triangles or tets) to connect the equispaced data points provided by PhysInterpToSimplexEquiSpaced. More...
 
void EquiSpacedToCoeffs (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs a projection/interpolation from the equispaced points sometimes used in post-processing onto the coefficient space. More...
 
template<class T >
boost::shared_ptr< T > as ()
 
void IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
- 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::GeomFactorsSharedPtr
v_GetMetricInfo () const
 
DNekMatSharedPtr BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
 
DNekMatSharedPtr BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd)
 
void ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs, std::vector< LibUtilities::BasisType > &fromType)
 
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 void v_NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, 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 void v_GetTracePhysVals (const int edge, const StdRegions::StdExpansionSharedPtr &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient)
 
virtual
StdRegions::StdExpansionSharedPtr 
v_GetStdExp (void) const
 
virtual
StdRegions::StdExpansionSharedPtr 
v_GetLinStdExp (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::BasisSharedPtr
v_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, std::vector< LibUtilities::BasisType > &fromType)
 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...
 
virtual void v_GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true)
 
virtual void v_ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion
DNekMatSharedPtr CreateStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr CreateStdStaticCondMatrix (const StdMatrixKey &mkey)
 Create the static condensation of a matrix when using a boundary interior decomposition. More...
 
IndexMapValuesSharedPtr CreateIndexMap (const IndexMapKey &ikey)
 Create an IndexMap which contains mapping information linking any specific element shape with either its boundaries, edges, faces, verteces, etc. More...
 
void BwdTrans_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GeneralMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
 
void LaplacianMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp_MatFree (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void HelmholtzMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
- 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::opLess
m_matrixManager
 
LibUtilities::NekManager
< MatrixKey, DNekScalBlkMat,
MatrixKey::opLess
m_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::BasisSharedPtr
m_base
 
int m_elmt_id
 
int m_ncoeffs
 
LibUtilities::NekManager
< StdMatrixKey, DNekMat,
StdMatrixKey::opLess
m_stdMatrixManager
 
LibUtilities::NekManager
< StdMatrixKey, DNekBlkMat,
StdMatrixKey::opLess
m_stdStaticCondMatrixManager
 
LibUtilities::NekManager
< IndexMapKey, IndexMapValues,
IndexMapKey::opLess
m_IndexMapManager
 
- 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 58 of file SegExp.cpp.

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

Copy Constructor.

Parameters
SExisting segment to duplicate.

Definition at line 79 of file SegExp.cpp.

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

Definition at line 94 of file SegExp.cpp.

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

Member Function Documentation

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

Definition at line 1208 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.

1209  {
1210  DNekScalMatSharedPtr returnval;
1211  NekDouble fac;
1213 
1214  ASSERTL2(m_metricinfo->GetGtype() !=
1216  "Geometric information is not set up");
1217 
1218  switch (mkey.GetMatrixType())
1219  {
1220  case StdRegions::eMass:
1221  {
1222  if ((m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1223  || (mkey.GetNVarCoeff()))
1224  {
1225  fac = 1.0;
1226  goto UseLocRegionsMatrix;
1227  }
1228  else
1229  {
1230  fac = (m_metricinfo->GetJac(ptsKeys))[0];
1231  goto UseStdRegionsMatrix;
1232  }
1233  }
1234  break;
1235  case StdRegions::eInvMass:
1236  {
1237  if ((m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1238  || (mkey.GetNVarCoeff()))
1239  {
1240  NekDouble one = 1.0;
1241  StdRegions::StdMatrixKey masskey(
1242  StdRegions::eMass,DetShapeType(), *this);
1243  DNekMatSharedPtr mat = GenMatrix(masskey);
1244  mat->Invert();
1245 
1246  returnval = MemoryManager<DNekScalMat>::
1247  AllocateSharedPtr(one,mat);
1248  }
1249  else
1250  {
1251  fac = 1.0/(m_metricinfo->GetJac(ptsKeys))[0];
1252  goto UseStdRegionsMatrix;
1253  }
1254  }
1255  break;
1259  {
1260  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed ||
1261  mkey.GetNVarCoeff())
1262  {
1263  fac = 1.0;
1264  goto UseLocRegionsMatrix;
1265  }
1266  else
1267  {
1268  int dir = 0;
1269  switch(mkey.GetMatrixType())
1270  {
1272  dir = 0;
1273  break;
1275  ASSERTL1(m_geom->GetCoordim() >= 2,
1276  "Cannot call eWeakDeriv2 in a "
1277  "coordinate system which is not at "
1278  "least two-dimensional");
1279  dir = 1;
1280  break;
1282  ASSERTL1(m_geom->GetCoordim() == 3,
1283  "Cannot call eWeakDeriv2 in a "
1284  "coordinate system which is not "
1285  "three-dimensional");
1286  dir = 2;
1287  break;
1288  default:
1289  break;
1290  }
1291 
1292  MatrixKey deriv0key(StdRegions::eWeakDeriv0,
1293  mkey.GetShapeType(), *this);
1294 
1295  DNekMatSharedPtr WeakDerivStd = GetStdMatrix(deriv0key);
1296  fac = m_metricinfo->GetDerivFactors(ptsKeys)[dir][0]*
1297  m_metricinfo->GetJac(ptsKeys)[0];
1298 
1299  returnval = MemoryManager<DNekScalMat>::
1300  AllocateSharedPtr(fac,WeakDerivStd);
1301  }
1302  }
1303  break;
1305  {
1306  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1307  {
1308  fac = 1.0;
1309  goto UseLocRegionsMatrix;
1310  }
1311  else
1312  {
1313  int coordim = m_geom->GetCoordim();
1314  fac = 0.0;
1315  for (int i = 0; i < coordim; ++i)
1316  {
1317  fac += m_metricinfo->GetDerivFactors(ptsKeys)[i][0]*
1318  m_metricinfo->GetDerivFactors(ptsKeys)[i][0];
1319  }
1320  fac *= m_metricinfo->GetJac(ptsKeys)[0];
1321  goto UseStdRegionsMatrix;
1322  }
1323  }
1324  break;
1326  {
1327  NekDouble factor =
1328  mkey.GetConstFactor(StdRegions::eFactorLambda);
1329  MatrixKey masskey(StdRegions::eMass,
1330  mkey.GetShapeType(), *this);
1331  DNekScalMat &MassMat = *(this->m_matrixManager[masskey]);
1332  MatrixKey lapkey(StdRegions::eLaplacian, mkey.GetShapeType(),
1333  *this, mkey.GetConstFactors(),
1334  mkey.GetVarCoeffs());
1335  DNekScalMat &LapMat = *(this->m_matrixManager[lapkey]);
1336 
1337  int rows = LapMat.GetRows();
1338  int cols = LapMat.GetColumns();
1339 
1340  DNekMatSharedPtr helm =
1342 
1343  NekDouble one = 1.0;
1344  (*helm) = LapMat + factor*MassMat;
1345 
1346  returnval =
1348  }
1349  break;
1354  {
1355  NekDouble one = 1.0;
1356 
1357  DNekMatSharedPtr mat = GenMatrix(mkey);
1358  returnval =
1360  }
1361  break;
1363  {
1364  NekDouble one = 1.0;
1365 
1366 // StdRegions::StdMatrixKey hkey(StdRegions::eHybridDGHelmholtz,
1367 // DetShapeType(),*this,
1368 // mkey.GetConstant(0),
1369 // mkey.GetConstant(1));
1370  MatrixKey hkey(StdRegions::eHybridDGHelmholtz,
1371  DetShapeType(),
1372  *this, mkey.GetConstFactors(),
1373  mkey.GetVarCoeffs());
1374  DNekMatSharedPtr mat = GenMatrix(hkey);
1375 
1376  mat->Invert();
1377  returnval =
1379  }
1380  break;
1382  {
1383  DNekMatSharedPtr m_Ix;
1384  Array<OneD, NekDouble> coords(1, 0.0);
1385  StdRegions::ConstFactorMap factors = mkey.GetConstFactors();
1386  int vertex = (int)factors[StdRegions::eFactorGaussVertex];
1387 
1388  coords[0] = (vertex == 0) ? -1.0 : 1.0;
1389 
1390  m_Ix = m_base[0]->GetI(coords);
1391  returnval =
1393  }
1394  break;
1395 
1396  UseLocRegionsMatrix:
1397  {
1398  DNekMatSharedPtr mat = GenMatrix(mkey);
1399  returnval =
1401  }
1402  break;
1403  UseStdRegionsMatrix:
1404  {
1405  DNekMatSharedPtr mat = GetStdMatrix(mkey);
1406  returnval =
1408  }
1409  break;
1410  default:
1411  NEKERROR(ErrorUtil::efatal, "Matrix creation not defined");
1412  break;
1413  }
1414 
1415  return returnval;
1416  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:470
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:191
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:242
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:131
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:252
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:130
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:706
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:250
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
Array< OneD, LibUtilities::BasisSharedPtr > m_base
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:250
DNekScalBlkMatSharedPtr Nektar::LocalRegions::SegExp::CreateStaticCondMatrix ( const MatrixKey mkey)
protected
Todo:
Really need a constructor which takes Arrays

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

1444  {
1445  DNekScalBlkMatSharedPtr returnval;
1446 
1447  ASSERTL2(m_metricinfo->GetGtype() !=
1449  "Geometric information is not set up");
1450 
1451  // set up block matrix system
1452  int nbdry = NumBndryCoeffs();
1453  int nint = m_ncoeffs - nbdry;
1454  Array<OneD, unsigned int> exp_size(2);
1455  exp_size[0] = nbdry;
1456  exp_size[1] = nint;
1457 
1458  /// \todo Really need a constructor which takes Arrays
1460  AllocateSharedPtr(exp_size,exp_size);
1461  NekDouble factor = 1.0;
1462 
1463  switch (mkey.GetMatrixType())
1464  {
1466  case StdRegions::eHelmholtz: // special case since Helmholtz
1467  // not defined in StdRegions
1468 
1469  // use Deformed case for both regular and deformed geometries
1470  factor = 1.0;
1471  goto UseLocRegionsMatrix;
1472  break;
1473  default:
1474  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1475  {
1476  factor = 1.0;
1477  goto UseLocRegionsMatrix;
1478  }
1479  else
1480  {
1481  DNekScalMatSharedPtr mat = GetLocMatrix(mkey);
1482  factor = mat->Scale();
1483  goto UseStdRegionsMatrix;
1484  }
1485  break;
1486  UseStdRegionsMatrix:
1487  {
1488  NekDouble invfactor = 1.0/factor;
1489  NekDouble one = 1.0;
1491  DNekScalMatSharedPtr Atmp;
1492  DNekMatSharedPtr Asubmat;
1493 
1494  returnval->SetBlock(0,0,Atmp =
1495  MemoryManager<DNekScalMat>::AllocateSharedPtr(
1496  factor,Asubmat = mat->GetBlock(0,0)));
1497  returnval->SetBlock(0,1,Atmp =
1498  MemoryManager<DNekScalMat>::AllocateSharedPtr(
1499  one,Asubmat = mat->GetBlock(0,1)));
1500  returnval->SetBlock(1,0,Atmp =
1501  MemoryManager<DNekScalMat>::AllocateSharedPtr(
1502  factor,Asubmat = mat->GetBlock(1,0)));
1503  returnval->SetBlock(1,1,Atmp =
1504  MemoryManager<DNekScalMat>::AllocateSharedPtr(
1505  invfactor,Asubmat = mat->GetBlock(1,1)));
1506  }
1507  break;
1508  UseLocRegionsMatrix:
1509  {
1510  int i,j;
1511  NekDouble invfactor = 1.0/factor;
1512  NekDouble one = 1.0;
1513  DNekScalMat &mat = *GetLocMatrix(mkey);
1516  DNekMatSharedPtr B =
1518  DNekMatSharedPtr C =
1520  DNekMatSharedPtr D =
1522 
1523  Array<OneD,unsigned int> bmap(nbdry);
1524  Array<OneD,unsigned int> imap(nint);
1525  GetBoundaryMap(bmap);
1526  GetInteriorMap(imap);
1527 
1528  for (i = 0; i < nbdry; ++i)
1529  {
1530  for (j = 0; j < nbdry; ++j)
1531  {
1532  (*A)(i,j) = mat(bmap[i],bmap[j]);
1533  }
1534 
1535  for (j = 0; j < nint; ++j)
1536  {
1537  (*B)(i,j) = mat(bmap[i],imap[j]);
1538  }
1539  }
1540 
1541  for (i = 0; i < nint; ++i)
1542  {
1543  for (j = 0; j < nbdry; ++j)
1544  {
1545  (*C)(i,j) = mat(imap[i],bmap[j]);
1546  }
1547 
1548  for (j = 0; j < nint; ++j)
1549  {
1550  (*D)(i,j) = mat(imap[i],imap[j]);
1551  }
1552  }
1553 
1554  // Calculate static condensed system
1555  if (nint)
1556  {
1557  D->Invert();
1558  (*B) = (*B)*(*D);
1559  (*A) = (*A) - (*B)*(*C);
1560  }
1561 
1562  DNekScalMatSharedPtr Atmp;
1563 
1564  returnval->SetBlock(0,0,Atmp = MemoryManager<DNekScalMat>::
1565  AllocateSharedPtr(factor,A));
1566  returnval->SetBlock(0,1,Atmp = MemoryManager<DNekScalMat>::
1567  AllocateSharedPtr(one,B));
1568  returnval->SetBlock(1,0,Atmp = MemoryManager<DNekScalMat>::
1569  AllocateSharedPtr(factor,C));
1570  returnval->SetBlock(1,1,Atmp = MemoryManager<DNekScalMat>::
1571  AllocateSharedPtr(invfactor,D));
1572  }
1573  }
1574 
1575 
1576  return returnval;
1577  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:131
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
DNekBlkMatSharedPtr GetStdStaticCondMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:711
boost::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:74
void GetInteriorMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:819
double NekDouble
boost::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:72
DNekScalMatSharedPtr GetLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:85
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:250
Geometry is curved or has non-constant factors.
NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag > DNekScalMat
void GetBoundaryMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:814
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 1631 of file SegExp.cpp.

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

1634  {
1635  // get Mass matrix inverse
1636  MatrixKey masskey(StdRegions::eInvMass,
1637  DetShapeType(),*this);
1638  DNekScalMatSharedPtr matsys = m_matrixManager[masskey];
1639 
1640  NekVector<NekDouble> in(m_ncoeffs,inarray,eCopy);
1641  NekVector<NekDouble> out(m_ncoeffs,outarray,eWrapper);
1642 
1643  out = (*matsys)*in;
1644  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:470
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: SegExp.h:250
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 1589 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().

1592  {
1593 
1594  int m;
1595  NekDouble sgn = 1;
1596 
1597  ASSERTL1(&inarray[0] != &outarray[0],
1598  "inarray and outarray can not be the same");
1599  switch(GetBasisType(0))
1600  {
1602  //Swap vertices
1603  outarray[0] = inarray[1];
1604  outarray[1] = inarray[0];
1605  // negate odd modes
1606  for(m = 2; m < m_ncoeffs; ++m)
1607  {
1608  outarray[m] = sgn*inarray[m];
1609  sgn = -sgn;
1610  }
1611  break;
1614  for(m = 0; m < m_ncoeffs; ++m)
1615  {
1616  outarray[m_ncoeffs-1-m] = inarray[m];
1617  }
1618  break;
1619  default:
1620  ASSERTL0(false,"This basis is not allowed in this method");
1621  break;
1622  }
1623  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
Principle Modified Functions .
Definition: BasisType.h:49
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:54
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:228
void Nektar::LocalRegions::SegExp::v_ComputeVertexNormal ( const int  vertex)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 895 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::StdExpansion1D::m_vertexNormals, and Vmath::Smul().

896  {
897  int i;
898  const SpatialDomains::GeomFactorsSharedPtr &geomFactors =
899  GetGeom()->GetMetricInfo();
900  SpatialDomains::GeomType type = geomFactors->GetGtype();
901  const Array<TwoD, const NekDouble> &gmat =
902  geomFactors->GetDerivFactors(GetPointsKeys());
903  int nqe = 1;
904  int vCoordDim = GetCoordim();
905 
907  Array<OneD, Array<OneD, NekDouble> >(vCoordDim);
908  Array<OneD, Array<OneD, NekDouble> > &normal =
910  for (i = 0; i < vCoordDim; ++i)
911  {
912  normal[i] = Array<OneD, NekDouble>(nqe);
913  }
914 
915  // Regular geometry case
916  if ((type == SpatialDomains::eRegular) ||
918  {
919  NekDouble vert;
920  // Set up normals
921  switch (vertex)
922  {
923  case 0:
924  for(i = 0; i < vCoordDim; ++i)
925  {
926  Vmath::Fill(nqe, -gmat[i][0], normal[i], 1);
927  }
928  break;
929  case 1:
930  for(i = 0; i < vCoordDim; ++i)
931  {
932  Vmath::Fill(nqe, gmat[i][0], normal[i], 1);
933  }
934  break;
935  default:
936  ASSERTL0(false,
937  "point is out of range (point < 2)");
938  }
939 
940  // normalise
941  vert = 0.0;
942  for (i =0 ; i < vCoordDim; ++i)
943  {
944  vert += normal[i][0]*normal[i][0];
945  }
946  vert = 1.0/sqrt(vert);
947  for (i = 0; i < vCoordDim; ++i)
948  {
949  Vmath::Smul(nqe, vert, normal[i], 1, normal[i], 1);
950  }
951  }
952  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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:213
double NekDouble
SpatialDomains::GeometrySharedPtr GetGeom() const
Definition: Expansion.cpp:161
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.
DNekMatSharedPtr Nektar::LocalRegions::SegExp::v_CreateStdMatrix ( const StdRegions::StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 1197 of file SegExp.cpp.

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

1199  {
1200  LibUtilities::BasisKey bkey = m_base[0]->GetBasisKey();
1203 
1204  return tmp->GetStdMatrix(mkey);
1205  }
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 1186 of file SegExp.cpp.

References m_staticCondMatrixManager.

1187  {
1188  m_staticCondMatrixManager.DeleteObject(mkey);
1189  }
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: SegExp.h:252
void Nektar::LocalRegions::SegExp::v_ExtractDataToCoeffs ( const NekDouble data,
const std::vector< unsigned int > &  nummodes,
const int  mode_offset,
NekDouble coeffs,
std::vector< LibUtilities::BasisType > &  fromType 
)
protectedvirtual

Unpack data from input file assuming it comes from.

Reimplemented from Nektar::LocalRegions::Expansion.

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

856  {
857  switch(m_base[0]->GetBasisType())
858  {
860  {
861  int fillorder = min((int) nummodes[mode_offset],m_ncoeffs);
862 
863  Vmath::Zero(m_ncoeffs,coeffs,1);
864  Vmath::Vcopy(fillorder,&data[0],1,&coeffs[0],1);
865  }
866  break;
868  {
869  // Assume that input is also Gll_Lagrange
870  // but no way to check;
871  LibUtilities::PointsKey p0(
872  nummodes[mode_offset],
875  p0,data, m_base[0]->GetPointsKey(), coeffs);
876  }
877  break;
879  {
880  // Assume that input is also Gauss_Lagrange
881  // but no way to check;
882  LibUtilities::PointsKey p0(
883  nummodes[mode_offset],
886  p0,data, m_base[0]->GetPointsKey(), coeffs);
887  }
888  break;
889  default:
890  ASSERTL0(false,
891  "basis is either not set up or not hierarchicial");
892  }
893  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
Principle Modified Functions .
Definition: BasisType.h:49
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:54
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:49
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
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:373
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:52
void Nektar::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 376 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().

379  {
380  if (m_base[0]->Collocation())
381  {
382  Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
383  }
384  else
385  {
386  v_IProductWRTBase(inarray,outarray);
387 
388  // get Mass matrix inverse
389  MatrixKey masskey(StdRegions::eInvMass, DetShapeType(), *this);
390  DNekScalMatSharedPtr matsys = m_matrixManager[masskey];
391 
392  // copy inarray in case inarray == outarray
393  NekVector<NekDouble> in(m_ncoeffs,outarray,eCopy);
394  NekVector<NekDouble> out(m_ncoeffs,outarray,eWrapper);
395 
396  out = (*matsys)*in;
397  }
398  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:470
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Inner product of inarray over region with respect to the expansion basis (this)->_Base[0] and return ...
Definition: SegExp.cpp:503
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:1061
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: SegExp.h:250
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 400 of file SegExp.cpp.

References ASSERTL0, ASSERTL1, Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::LibUtilities::eGauss_Lagrange, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGLL_Lagrange, Nektar::StdRegions::eMass, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::GetPointsType(), 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().

403  {
404  if(m_base[0]->Collocation())
405  {
406  Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
407  }
408  else
409  {
410  int nInteriorDofs = m_ncoeffs-2;
411  int offset;
412 
413  switch (m_base[0]->GetBasisType())
414  {
416  {
417  offset = 1;
418  }
419  break;
421  {
422  nInteriorDofs = m_ncoeffs;
423  offset = 0;
424  }
425  break;
428  {
429  ASSERTL1(m_base[0]->GetPointsType() == LibUtilities::eGaussLobattoLegendre,"Cannot use FwdTrans_BndConstrained method with non GLL points");
430  offset = 2;
431  }
432  break;
433  default:
434  ASSERTL0(false,"This type of FwdTrans is not defined"
435  "for this expansion type");
436  }
437 
438  fill(outarray.get(), outarray.get()+m_ncoeffs, 0.0 );
439 
441  {
442 
443  outarray[GetVertexMap(0)] = inarray[0];
444  outarray[GetVertexMap(1)] =
445  inarray[m_base[0]->GetNumPoints()-1];
446 
447  if (m_ncoeffs>2)
448  {
449  // ideally, we would like to have tmp0 to be replaced
450  // by outarray (currently MassMatrixOp does not allow
451  // aliasing)
452  Array<OneD, NekDouble> tmp0(m_ncoeffs);
453  Array<OneD, NekDouble> tmp1(m_ncoeffs);
454 
455  StdRegions::StdMatrixKey stdmasskey(
457  MassMatrixOp(outarray,tmp0,stdmasskey);
458  v_IProductWRTBase(inarray,tmp1);
459 
460  Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
461 
462  // get Mass matrix inverse (only of interior DOF)
463  MatrixKey masskey(
465  DNekScalMatSharedPtr matsys =
466  (m_staticCondMatrixManager[masskey])->GetBlock(1,1);
467 
468  Blas::Dgemv('N',nInteriorDofs,nInteriorDofs,
469  matsys->Scale(),
470  &((matsys->GetOwnedMatrix())->GetPtr())[0],
471  nInteriorDofs,tmp1.get()+offset,1,0.0,
472  outarray.get()+offset,1);
473  }
474  }
475  else
476  {
477  SegExp::v_FwdTrans(inarray, outarray);
478 
479  }
480  }
481  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:470
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:976
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: SegExp.h:252
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:376
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:503
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
int GetVertexMap(const int localVertexId, bool useCoeffPacking=false)
Definition: StdExpansion.h:824
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:343
Lagrange for SEM basis .
Definition: BasisType.h:53
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:216
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:52
DNekMatSharedPtr Nektar::LocalRegions::SegExp::v_GenMatrix ( const StdRegions::StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

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

1420  {
1421  DNekMatSharedPtr returnval;
1422 
1423  switch (mkey.GetMatrixType())
1424  {
1431  returnval = Expansion1D::v_GenMatrix(mkey);
1432  break;
1433  default:
1434  returnval = StdSegExp::v_GenMatrix(mkey);
1435  break;
1436  }
1437 
1438  return returnval;
1439  }
virtual DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey)
Definition: Expansion1D.cpp:45
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
const LibUtilities::BasisSharedPtr & Nektar::LocalRegions::SegExp::v_GetBasis ( int  dir) const
protectedvirtual

Definition at line 830 of file SegExp.cpp.

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

831  {
832  return GetBasis(dir);
833  }
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 679 of file SegExp.cpp.

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

682  {
683  int i;
684 
685  ASSERTL1(Lcoords[0] >= -1.0&& Lcoords[0] <= 1.0,
686  "Local coordinates are not in region [-1,1]");
687 
688  m_geom->FillGeom();
689  for(i = 0; i < m_geom->GetCoordim(); ++i)
690  {
691  coords[i] = m_geom->GetCoord(i,Lcoords);
692  }
693  }
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:130
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
int Nektar::LocalRegions::SegExp::v_GetCoordim ( void  )
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion1D.

Definition at line 804 of file SegExp.cpp.

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

805  {
806  return m_geom->GetCoordim();
807  }
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:130
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 695 of file SegExp.cpp.

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

699  {
700  Expansion::v_GetCoords(coords_0, coords_1, coords_2);
701  }
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3)
Definition: Expansion.cpp:224
StdRegions::StdExpansionSharedPtr Nektar::LocalRegions::SegExp::v_GetLinStdExp ( void  ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 795 of file SegExp.cpp.

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

796  {
797  LibUtilities::BasisKey bkey0(m_base[0]->GetBasisType(),
798  2, m_base[0]->GetPointsKey());
799 
801  ::AllocateSharedPtr( bkey0);
802  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
Array< OneD, LibUtilities::BasisSharedPtr > m_base
DNekScalMatSharedPtr Nektar::LocalRegions::SegExp::v_GetLocMatrix ( const MatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1191 of file SegExp.cpp.

References m_matrixManager.

1192  {
1193  return m_matrixManager[mkey];
1194  }
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: SegExp.h:250
DNekScalBlkMatSharedPtr Nektar::LocalRegions::SegExp::v_GetLocStaticCondMatrix ( const MatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1180 of file SegExp.cpp.

References m_staticCondMatrixManager.

1182  {
1183  return m_staticCondMatrixManager[mkey];
1184  }
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: SegExp.h:252
int Nektar::LocalRegions::SegExp::v_GetNcoeffs ( void  ) const
protectedvirtual

Definition at line 825 of file SegExp.cpp.

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

826  {
827  return m_ncoeffs;
828  }
int Nektar::LocalRegions::SegExp::v_GetNumPoints ( const int  dir) const
protectedvirtual

Definition at line 820 of file SegExp.cpp.

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

821  {
822  return GetNumPoints(dir);
823  }
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 809 of file SegExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 783 of file SegExp.cpp.

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

784  {
785  return m_geom->GetPorient(point);
786  }
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:130
bg::model::point< double, 3, bg::cs::cartesian > point
Definition: BLMesh.cpp:53
StdRegions::StdExpansionSharedPtr Nektar::LocalRegions::SegExp::v_GetStdExp ( void  ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 789 of file SegExp.cpp.

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

790  {
792  ::AllocateSharedPtr(m_base[0]->GetBasisKey());
793  }
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_GetTracePhysVals ( const int  edge,
const StdRegions::StdExpansionSharedPtr EdgeExp,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
StdRegions::Orientation  orient 
)
protectedvirtual

Definition at line 740 of file SegExp.cpp.

References v_GetVertexPhysVals().

746  {
747  NekDouble result;
748  v_GetVertexPhysVals(edge, inarray, result);
749  outarray[0] = result;
750  }
double NekDouble
virtual void v_GetVertexPhysVals(const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray)
Definition: SegExp.cpp:704
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 704 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.

Referenced by v_GetTracePhysVals().

708  {
709  int nquad = m_base[0]->GetNumPoints();
710 
712  {
713  switch (vertex)
714  {
715  case 0:
716  outarray = inarray[0];
717  break;
718  case 1:
719  outarray = inarray[nquad - 1];
720  break;
721  }
722  }
723  else
724  {
727 
728  StdRegions::StdMatrixKey key(
730  DetShapeType(),*this,factors);
731 
732  DNekScalMatSharedPtr mat_gauss = m_matrixManager[key];
733 
734  outarray = Blas::Ddot(nquad, mat_gauss->GetOwnedMatrix()
735  ->GetPtr().get(), 1, &inarray[0], 1);
736  }
737  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:470
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:252
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:49
T Ddot(int n, const Array< OneD, const T > &w, const int incw, const Array< OneD, const T > &x, const int incx, const Array< OneD, const int > &y, const int incy)
Definition: VmathArray.hpp:436
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:250
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 1065 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().

1069  {
1070  int nquad = m_base[0]->GetNumPoints();
1071  const Array<TwoD, const NekDouble>& gmat =
1072  m_metricinfo->GetDerivFactors(GetPointsKeys());
1073  const NekDouble lambda =
1074  mkey.GetConstFactor(StdRegions::eFactorLambda);
1075 
1076  Array<OneD,NekDouble> physValues(nquad);
1077  Array<OneD,NekDouble> dPhysValuesdx(nquad);
1078  Array<OneD,NekDouble> wsp(m_ncoeffs);
1079 
1080  BwdTrans(inarray, physValues);
1081 
1082  // mass matrix operation
1083  v_IProductWRTBase((m_base[0]->GetBdata()),physValues,wsp,1);
1084 
1085  // Laplacian matrix operation
1086  switch (m_geom->GetCoordim())
1087  {
1088  case 1:
1089  {
1090  PhysDeriv(physValues,dPhysValuesdx);
1091 
1092  // multiply with the proper geometric factors
1093  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1094  {
1095  Vmath::Vmul(nquad,
1096  &gmat[0][0],1,dPhysValuesdx.get(),1,
1097  dPhysValuesdx.get(),1);
1098  }
1099  else
1100  {
1101  Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.get(), 1);
1102  }
1103  }
1104  break;
1105  case 2:
1106  {
1107  Array<OneD,NekDouble> dPhysValuesdy(nquad);
1108 
1109  PhysDeriv(physValues, dPhysValuesdx, dPhysValuesdy);
1110 
1111  // multiply with the proper geometric factors
1112  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1113  {
1114  Vmath::Vmul (nquad,
1115  &gmat[0][0], 1, dPhysValuesdx.get(), 1,
1116  dPhysValuesdx.get(), 1);
1117  Vmath::Vvtvp(nquad,
1118  &gmat[1][0], 1, dPhysValuesdy.get(), 1,
1119  dPhysValuesdx.get(), 1,
1120  dPhysValuesdx.get(), 1);
1121  }
1122  else
1123  {
1124  Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.get(), 1);
1125  Blas::Daxpy(nquad,
1126  gmat[1][0], dPhysValuesdy.get(), 1,
1127  dPhysValuesdx.get(), 1);
1128  }
1129  }
1130  break;
1131  case 3:
1132  {
1133  Array<OneD,NekDouble> dPhysValuesdy(nquad);
1134  Array<OneD,NekDouble> dPhysValuesdz(nquad);
1135 
1136  PhysDeriv(physValues, dPhysValuesdx,
1137  dPhysValuesdy, dPhysValuesdz);
1138 
1139  // multiply with the proper geometric factors
1140  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1141  {
1142  Vmath::Vmul (nquad,
1143  &gmat[0][0], 1, dPhysValuesdx.get(), 1,
1144  dPhysValuesdx.get(), 1);
1145  Vmath::Vvtvp(nquad,
1146  &gmat[1][0], 1, dPhysValuesdy.get(), 1,
1147  dPhysValuesdx.get(), 1,
1148  dPhysValuesdx.get(), 1);
1149  Vmath::Vvtvp(nquad,
1150  &gmat[2][0], 1, dPhysValuesdz.get(), 1,
1151  dPhysValuesdx.get(), 1,
1152  dPhysValuesdx.get(), 1);
1153  }
1154  else
1155  {
1156  Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.get(), 1);
1157  Blas::Daxpy(nquad,
1158  gmat[1][0], dPhysValuesdy.get(), 1,
1159  dPhysValuesdx.get(), 1);
1160  Blas::Daxpy(nquad,
1161  gmat[2][0], dPhysValuesdz.get(),
1162  1, dPhysValuesdx.get(), 1);
1163  }
1164  }
1165  break;
1166  default:
1167  ASSERTL0(false,"Wrong number of dimensions");
1168  break;
1169  }
1170 
1171  v_IProductWRTBase(m_base[0]->GetDbdata(),dPhysValuesdx,outarray,1);
1172  Blas::Daxpy(m_ncoeffs, lambda, wsp.get(), 1, outarray.get(), 1);
1173  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:442
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:131
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:130
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:503
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:531
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:183
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 117 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().

119  {
120  int nquad0 = m_base[0]->GetNumPoints();
121  Array<OneD, const NekDouble> jac = m_metricinfo->GetJac(GetPointsKeys());
122  NekDouble ival;
123  Array<OneD,NekDouble> tmp(nquad0);
124 
125  // multiply inarray with Jacobian
126  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
127  {
128  Vmath::Vmul(nquad0, jac, 1, inarray, 1, tmp,1);
129  }
130  else
131  {
132  Vmath::Smul(nquad0, jac[0], inarray, 1, tmp, 1);
133  }
134 
135  // call StdSegExp version;
136  ival = StdSegExp::v_Integral(tmp);
137  //ival = StdSegExp::Integral(tmp);
138  return ival;
139  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:131
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:213
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:183
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 503 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().

506  {
507  v_IProductWRTBase(m_base[0]->GetBdata(),inarray,outarray,1);
508  }
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:503
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 538 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().

543  {
544  int nquad0 = m_base[0]->GetNumPoints();
545  Array<OneD, const NekDouble> jac = m_metricinfo->GetJac(GetPointsKeys());
546  Array<OneD,NekDouble> tmp(nquad0);
547 
548 
549  // multiply inarray with Jacobian
550  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
551  {
552  Vmath::Vmul(nquad0, jac, 1, inarray, 1, tmp, 1);
553  }
554  else
555  {
556  Vmath::Smul(nquad0, jac[0], inarray, 1, tmp, 1);
557  }
558  StdSegExp::v_IProductWRTBase(base,tmp,outarray,coll_check);
559  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:131
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:213
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:183
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 562 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().

566  {
567  int nquad = m_base[0]->GetNumPoints();
568  const Array<TwoD, const NekDouble>& gmat =
569  m_metricinfo->GetDerivFactors(GetPointsKeys());
570 
571  Array<OneD, NekDouble> tmp1(nquad);
572 
573  switch(dir)
574  {
575  case 0:
576  {
577  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
578  {
579  Vmath::Vmul(nquad,gmat[0],1,inarray,1,tmp1,1);
580  }
581  else
582  {
583  Vmath::Smul(nquad, gmat[0][0], inarray, 1, tmp1, 1);
584  }
585  }
586  break;
587  case 1:
588  {
589  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
590  {
591  Vmath::Vmul(nquad,gmat[1],1,inarray,1,tmp1,1);
592  }
593  else
594  {
595  Vmath::Smul(nquad, gmat[1][0], inarray, 1, tmp1, 1);
596  }
597  }
598  break;
599  case 2:
600  {
601  ASSERTL1(m_geom->GetCoordim() == 3,"input dir is out of range");
602  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
603  {
604  Vmath::Vmul(nquad,gmat[2],1,inarray,1,tmp1,1);
605  }
606  else
607  {
608  Vmath::Smul(nquad, gmat[2][0], inarray, 1, tmp1, 1);
609  }
610  }
611  break;
612  default:
613  {
614  ASSERTL1(false,"input dir is out of range");
615  }
616  break;
617  }
618  v_IProductWRTBase(m_base[0]->GetDbdata(),tmp1,outarray,1);
619  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:131
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:130
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:503
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:213
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
Array< OneD, LibUtilities::BasisSharedPtr > m_base
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:183
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 959 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().

963  {
964  int nquad = m_base[0]->GetNumPoints();
965  const Array<TwoD, const NekDouble>& gmat =
966  m_metricinfo->GetDerivFactors(GetPointsKeys());
967 
968  Array<OneD,NekDouble> physValues(nquad);
969  Array<OneD,NekDouble> dPhysValuesdx(nquad);
970 
971  BwdTrans(inarray,physValues);
972 
973  // Laplacian matrix operation
974  switch (m_geom->GetCoordim())
975  {
976  case 1:
977  {
978  PhysDeriv(physValues,dPhysValuesdx);
979 
980  // multiply with the proper geometric factors
981  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
982  {
983  Vmath::Vmul(nquad,
984  &gmat[0][0],1,dPhysValuesdx.get(),1,
985  dPhysValuesdx.get(),1);
986  }
987  else
988  {
989  Blas::Dscal(nquad,
990  gmat[0][0], dPhysValuesdx.get(), 1);
991  }
992  }
993  break;
994  case 2:
995  {
996  Array<OneD,NekDouble> dPhysValuesdy(nquad);
997 
998  PhysDeriv(physValues,dPhysValuesdx,dPhysValuesdy);
999 
1000  // multiply with the proper geometric factors
1001  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1002  {
1003  Vmath::Vmul (nquad,
1004  &gmat[0][0],1,dPhysValuesdx.get(),1,
1005  dPhysValuesdx.get(),1);
1006  Vmath::Vvtvp(nquad,
1007  &gmat[1][0],1,dPhysValuesdy.get(),1,
1008  dPhysValuesdx.get(),1,
1009  dPhysValuesdx.get(),1);
1010  }
1011  else
1012  {
1013  Blas::Dscal(nquad,
1014  gmat[0][0], dPhysValuesdx.get(), 1);
1015  Blas::Daxpy(nquad,
1016  gmat[1][0], dPhysValuesdy.get(), 1,
1017  dPhysValuesdx.get(), 1);
1018  }
1019  }
1020  break;
1021  case 3:
1022  {
1023  Array<OneD,NekDouble> dPhysValuesdy(nquad);
1024  Array<OneD,NekDouble> dPhysValuesdz(nquad);
1025 
1026  PhysDeriv(physValues,dPhysValuesdx,
1027  dPhysValuesdy,dPhysValuesdz);
1028 
1029  // multiply with the proper geometric factors
1030  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1031  {
1032  Vmath::Vmul (nquad,
1033  &gmat[0][0], 1, dPhysValuesdx.get(), 1,
1034  dPhysValuesdx.get(),1);
1035  Vmath::Vvtvp(nquad,
1036  &gmat[1][0], 1, dPhysValuesdy.get(), 1,
1037  dPhysValuesdx.get(),1,
1038  dPhysValuesdx.get(),1);
1039  Vmath::Vvtvp(nquad,
1040  &gmat[2][0], 1, dPhysValuesdz.get(), 1,
1041  dPhysValuesdx.get(),1,
1042  dPhysValuesdx.get(),1);
1043  }
1044  else
1045  {
1046  Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.get(), 1);
1047  Blas::Daxpy(nquad,
1048  gmat[1][0], dPhysValuesdy.get(), 1,
1049  dPhysValuesdx.get(), 1);
1050  Blas::Daxpy(nquad,
1051  gmat[2][0], dPhysValuesdz.get(), 1,
1052  dPhysValuesdx.get(), 1);
1053  }
1054  }
1055  break;
1056  default:
1057  ASSERTL0(false,"Wrong number of dimensions");
1058  break;
1059  }
1060 
1061  v_IProductWRTBase(m_base[0]->GetDbdata(),dPhysValuesdx,outarray,1);
1062  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:442
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:131
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:130
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:503
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:531
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:183
SpatialDomains::GeomType Nektar::LocalRegions::SegExp::v_MetricInfoType ( )
protectedvirtual

Definition at line 815 of file SegExp.cpp.

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

816  {
817  return m_metricinfo->GetGtype();
818  }
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:131
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 621 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().

625  {
626  int nq = m_base[0]->GetNumPoints();
627  Array<OneD, NekDouble > Fn(nq);
628 // cout << "I am segment " << GetGeom()->GetGlobalID() << endl;
629 // cout << "I want edge " << GetLeftAdjacentElementEdge() << endl;
630 // @TODO: This routine no longer makes sense as a normal is not unique to an edge
631  const Array<OneD, const Array<OneD, NekDouble> >
632  &normals =
635  Vmath::Vmul (nq, &Fx[0], 1, &normals[0][0], 1, &Fn[0], 1);
636  Vmath::Vvtvp(nq, &Fy[0], 1, &normals[1][0], 1, &Fn[0], 1, &Fn[0], 1);
637 
638  v_IProductWRTBase(Fn,outarray);
639  }
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:442
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:503
Expansion2DSharedPtr GetLeftAdjacentElementExp() const
Definition: Expansion1D.h:123
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:183
void Nektar::LocalRegions::SegExp::v_NormVectorIProductWRTBase ( const Array< OneD, const Array< OneD, NekDouble > > &  Fvec,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 641 of file SegExp.cpp.

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

644  {
645  NormVectorIProductWRTBase(Fvec[0], Fvec[1], outarray);
646  }
void NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:733
int Nektar::LocalRegions::SegExp::v_NumBndryCoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 836 of file SegExp.cpp.

837  {
838  return 2;
839  }
int Nektar::LocalRegions::SegExp::v_NumDGBndryCoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 842 of file SegExp.cpp.

843  {
844  return 2;
845  }
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 165 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().

170  {
171  int nquad0 = m_base[0]->GetNumPoints();
172  Array<TwoD, const NekDouble> gmat =
173  m_metricinfo->GetDerivFactors(GetPointsKeys());
174  Array<OneD,NekDouble> diff(nquad0);
175 
176  //StdExpansion1D::PhysTensorDeriv(inarray,diff);
177  PhysTensorDeriv(inarray,diff);
178  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
179  {
180  if(out_d0.num_elements())
181  {
182  Vmath::Vmul(nquad0,&gmat[0][0],1,&diff[0],1,
183  &out_d0[0],1);
184  }
185 
186  if(out_d1.num_elements())
187  {
188  Vmath::Vmul(nquad0,&gmat[1][0],1,&diff[0],1,
189  &out_d1[0],1);
190  }
191 
192  if(out_d2.num_elements())
193  {
194  Vmath::Vmul(nquad0,&gmat[2][0],1,&diff[0],1,
195  &out_d2[0],1);
196  }
197  }
198  else
199  {
200  if(out_d0.num_elements())
201  {
202  Vmath::Smul(nquad0, gmat[0][0], diff, 1,
203  out_d0, 1);
204  }
205 
206  if(out_d1.num_elements())
207  {
208  Vmath::Smul(nquad0, gmat[1][0], diff, 1,
209  out_d1, 1);
210  }
211 
212  if(out_d2.num_elements())
213  {
214  Vmath::Smul(nquad0, gmat[2][0], diff, 1,
215  out_d2, 1);
216  }
217  }
218  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:131
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:213
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:183
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 318 of file SegExp.cpp.

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

321  {
322  switch(dir)
323  {
324  case 0:
325  {
326  PhysDeriv(inarray, outarray, NullNekDouble1DArray,
328  }
329  break;
330  case 1:
331  {
332  PhysDeriv(inarray, NullNekDouble1DArray, outarray,
334  }
335  break;
336  case 2:
337  {
339  NullNekDouble1DArray, outarray);
340  }
341  break;
342  default:
343  {
344  ASSERTL1(false,"input dir is out of range");
345  }
346  break;
347  }
348  }
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:228
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 268 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().

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

231  {
232  int nquad0 = m_base[0]->GetNumPoints();
233  int coordim = m_geom->GetCoordim();
234  Array<OneD, NekDouble> diff (nquad0);
235  //this operation is needed if you put out_ds==inarray
236  Vmath::Zero(nquad0,out_ds,1);
237  switch(coordim)
238  {
239  case 2:
240  //diff= dU/de
241  Array<OneD,NekDouble> diff(nquad0);
242 
243  PhysTensorDeriv(inarray,diff);
244 
245  //get dS/de= (Jac)^-1
246  Array<OneD, NekDouble> Jac = m_metricinfo->GetJac(GetPointsKeys());
247  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
248  {
249  //calculate the derivative as (dU/de)*(Jac)^-1
250  Vmath::Vdiv(nquad0,diff,1,Jac ,1,out_ds,1);
251  }
252  else
253  {
254  NekDouble invJac = 1/Jac[0];
255  Vmath::Smul(nquad0, invJac,diff,1,out_ds,1);
256  }
257  }
258  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:131
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:241
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:130
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:213
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:373
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 666 of file SegExp.cpp.

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

669  {
670  Array<OneD,NekDouble> Lcoord = Array<OneD,NekDouble>(1);
671 
672  ASSERTL0(m_geom,"m_geom not defined");
673  m_geom->GetLocCoords(coord,Lcoord);
674 
675  return StdSegExp::v_PhysEvaluate(Lcoord, physvals);
676  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:130
void Nektar::LocalRegions::SegExp::v_SetCoeffsToOrientation ( Array< OneD, NekDouble > &  coeffs,
StdRegions::Orientation  dir 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 756 of file SegExp.cpp.

759  {
760  v_SetCoeffsToOrientation(dir,coeffs,coeffs);
761  }
virtual void v_SetCoeffsToOrientation(Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
Definition: SegExp.cpp:756
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 763 of file SegExp.cpp.

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

767  {
768 
769  if (dir == StdRegions::eBackwards)
770  {
771  if (&inarray[0] != &outarray[0])
772  {
773  Array<OneD,NekDouble> intmp (inarray);
774  ReverseCoeffsAndSign(intmp,outarray);
775  }
776  else
777  {
778  ReverseCoeffsAndSign(inarray,outarray);
779  }
780  }
781  }
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:1589
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 658 of file SegExp.cpp.

661  {
662  // Evaluate point in local (eta) coordinates.
663  return StdSegExp::v_PhysEvaluate(Lcoord,physvals);
664  }

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