Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator 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
 
int GetShapeDimension () const
 
bool IsBoundaryInteriorExpansion ()
 
bool IsNodalNonTensorialExp ()
 
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space. More...
 
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space. More...
 
void FwdTrans_BndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble Integral (const Array< OneD, const NekDouble > &inarray)
 This function integrates the specified function over the domain. More...
 
void FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 This function fills the array outarray with the mode-th mode of the expansion. More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 this function calculates the inner product of a given function f with the different modes of the expansion More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &base, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, int coll_check)
 
void IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
int GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void SetElmtId (const int id)
 Set the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
void GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2=NullNekDouble1DArray, Array< OneD, NekDouble > &coords_3=NullNekDouble1DArray)
 this function returns the physical coordinates of the quadrature points of the expansion More...
 
void GetCoord (const Array< OneD, const NekDouble > &Lcoord, Array< OneD, NekDouble > &coord)
 given the coordinates of a point of the element in the local collapsed coordinate system, this function calculates the physical coordinates of the point More...
 
DNekMatSharedPtr GetStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr GetStdStaticCondMatrix (const StdMatrixKey &mkey)
 
IndexMapValuesSharedPtr GetIndexMap (const IndexMapKey &ikey)
 
const Array< OneD, const
NekDouble > & 
GetPhysNormals (void)
 
void SetPhysNormals (Array< OneD, const NekDouble > &normal)
 
virtual void SetUpPhysNormals (const int edge)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
 
void NormVectorIProductWRTBase (const Array< OneD, const Array< OneD, NekDouble > > &Fvec, Array< OneD, NekDouble > &outarray)
 
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
void DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
 
StdRegions::Orientation GetForient (int face)
 
StdRegions::Orientation GetEorient (int edge)
 
StdRegions::Orientation GetPorient (int point)
 
StdRegions::Orientation GetCartesianEorient (int edge)
 
void SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
 
void SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
int CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
void ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs)
 
NekDouble StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
 
int GetCoordim ()
 
void GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
void GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
int GetVertexMap (const int localVertexId, bool useCoeffPacking=false)
 
void GetEdgeInteriorMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
 
void GetFaceInteriorMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
 
void GetEdgeToElementMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int P=-1)
 
void GetFaceToElementMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int nummodesA=-1, int nummodesB=-1)
 
void GetEdgePhysVals (const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Extract the physical values along edge edge from inarray into outarray following the local edge orientation and point distribution defined by defined in EdgeExp. More...
 
void GetEdgePhysVals (const int edge, const boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetTracePhysVals (const int edge, const boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetVertexPhysVals (const int vertex, const Array< OneD, const NekDouble > &inarray, NekDouble &outarray)
 
void GetEdgeInterpVals (const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetEdgeQFactors (const int edge, Array< OneD, NekDouble > &outarray)
 Extract the metric factors to compute the contravariant fluxes along edge edge and stores them into outarray following the local edge orientation (i.e. anticlockwise convention). More...
 
void GetFacePhysVals (const int face, const boost::shared_ptr< StdExpansion > &FaceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=eNoOrientation)
 
void GetEdgePhysMap (const int edge, Array< OneD, int > &outarray)
 
void GetFacePhysMap (const int face, Array< OneD, int > &outarray)
 
void MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
DNekMatSharedPtr CreateGeneralMatrix (const StdMatrixKey &mkey)
 this function generates the mass matrix $\mathbf{M}[i][j] = \int \phi_i(\mathbf{x}) \phi_j(\mathbf{x}) d\mathbf{x}$ More...
 
void GeneralMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
DNekMatSharedPtr GenMatrix (const StdMatrixKey &mkey)
 
void PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void PhysDeriv_s (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_ds)
 
void PhysDeriv_n (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_dn)
 
void PhysDirectionalDeriv (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &outarray)
 
void StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
void StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void AddRobinMassMatrix (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
void AddRobinEdgeContribution (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, Array< OneD, NekDouble > &coeffs)
 
NekDouble PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
NekDouble PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain. More...
 
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta. More...
 
const boost::shared_ptr
< SpatialDomains::GeomFactors > & 
GetMetricInfo (void) const
 
virtual int v_GetElmtId ()
 Get the element id of this expansion when used in a list by returning value of m_elmt_id. More...
 
virtual 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 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 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)
 Unpack data from input file assuming it comes from. More...
 
virtual const Array< OneD,
const NekDouble > & 
v_GetPhysNormals (void)
 
virtual void v_LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
 
virtual void v_HelmholtzMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
 
virtual DNekMatSharedPtr v_GenMatrix (const StdRegions::StdMatrixKey &mkey)
 
DNekScalMatSharedPtr CreateMatrix (const MatrixKey &mkey)
 
virtual DNekMatSharedPtr v_CreateStdMatrix (const StdRegions::StdMatrixKey &mkey)
 
DNekScalBlkMatSharedPtr CreateStaticCondMatrix (const MatrixKey &mkey)
 
virtual DNekScalMatSharedPtr v_GetLocMatrix (const MatrixKey &mkey)
 
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const MatrixKey &mkey)
 
void v_DropLocStaticCondMatrix (const MatrixKey &mkey)
 
- Protected Member Functions inherited from Nektar::StdRegions::StdSegExp
virtual void v_StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 
virtual void v_StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Backward transform from coefficient space given in inarray and evaluate at the physical quadrature points outarray. More...
 
virtual void v_BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
 
virtual void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 
virtual void v_GetBoundaryMap (Array< OneD, unsigned int > &outarray)
 
virtual void v_GetInteriorMap (Array< OneD, unsigned int > &outarray)
 
virtual int v_GetVertexMap (int localVertexId, bool useCoeffPacking=false)
 
virtual int v_GetNverts () const
 
virtual bool v_IsBoundaryInteriorExpansion ()
 
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
 
virtual LibUtilities::ShapeType v_DetShapeType () const
 Return Shape of region, using ShapeType enum list. i.e. Segment. More...
 
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion
DNekMatSharedPtr CreateStdMatrix (const StdMatrixKey &mkey)
 
DNekBlkMatSharedPtr CreateStdStaticCondMatrix (const StdMatrixKey &mkey)
 Create the static condensation of a matrix when using a boundary interior decomposition. More...
 
IndexMapValuesSharedPtr CreateIndexMap (const IndexMapKey &ikey)
 Create an IndexMap which contains mapping information linking any specific element shape with either its boundaries, edges, faces, verteces, etc. More...
 
void BwdTrans_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GeneralMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
 
void LaplacianMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LaplacianMatrixOp_MatFree (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDerivMatrixOp_MatFree (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void WeakDirectionalDerivMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void MassLevelCurvatureMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void LinearAdvectionDiffusionReactionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
 
void HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
void HelmholtzMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
 
- Protected Member Functions inherited from Nektar::LocalRegions::Expansion1D
virtual void v_AddRobinMassMatrix (const int vert, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
 
virtual void v_AddRobinEdgeContribution (const int vert, const Array< OneD, const NekDouble > &primCoeffs, Array< OneD, NekDouble > &coeffs)
 
- Protected Member Functions inherited from Nektar::LocalRegions::Expansion
void ComputeLaplacianMetric ()
 
void ComputeQuadratureMetric ()
 
virtual void v_MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_ComputeLaplacianMetric ()
 
virtual DNekMatSharedPtr v_BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
 
virtual DNekMatSharedPtr v_BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd)
 
virtual void v_AddEdgeNormBoundaryInt (const int edge, const boost::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddEdgeNormBoundaryInt (const int edge, const boost::shared_ptr< Expansion > &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddFaceNormBoundaryInt (const int face, const boost::shared_ptr< Expansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
virtual void v_DGDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &coeffs, Array< OneD, NekDouble > &outarray)
 

Private Member Functions

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

Private Attributes

LibUtilities::NekManager
< MatrixKey, DNekScalMat,
MatrixKey::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 57 of file SegExp.cpp.

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

Copy Constructor.

Parameters
SExisting segment to duplicate.

Definition at line 78 of file SegExp.cpp.

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

Definition at line 93 of file SegExp.cpp.

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

Member Function Documentation

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

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

1198  {
1199  DNekScalMatSharedPtr returnval;
1200  NekDouble fac;
1202 
1203  ASSERTL2(m_metricinfo->GetGtype() !=
1205  "Geometric information is not set up");
1206 
1207  switch (mkey.GetMatrixType())
1208  {
1209  case StdRegions::eMass:
1210  {
1211  if ((m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1212  || (mkey.GetNVarCoeff()))
1213  {
1214  fac = 1.0;
1215  goto UseLocRegionsMatrix;
1216  }
1217  else
1218  {
1219  fac = (m_metricinfo->GetJac(ptsKeys))[0];
1220  goto UseStdRegionsMatrix;
1221  }
1222  }
1223  break;
1224  case StdRegions::eInvMass:
1225  {
1226  if ((m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1227  || (mkey.GetNVarCoeff()))
1228  {
1229  NekDouble one = 1.0;
1230  StdRegions::StdMatrixKey masskey(
1231  StdRegions::eMass,DetShapeType(), *this);
1232  DNekMatSharedPtr mat = GenMatrix(masskey);
1233  mat->Invert();
1234 
1235  returnval = MemoryManager<DNekScalMat>::
1236  AllocateSharedPtr(one,mat);
1237  }
1238  else
1239  {
1240  fac = 1.0/(m_metricinfo->GetJac(ptsKeys))[0];
1241  goto UseStdRegionsMatrix;
1242  }
1243  }
1244  break;
1248  {
1249  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed ||
1250  mkey.GetNVarCoeff())
1251  {
1252  fac = 1.0;
1253  goto UseLocRegionsMatrix;
1254  }
1255  else
1256  {
1257  int dir = 0;
1258  switch(mkey.GetMatrixType())
1259  {
1261  dir = 0;
1262  break;
1264  ASSERTL1(m_geom->GetCoordim() >= 2,
1265  "Cannot call eWeakDeriv2 in a "
1266  "coordinate system which is not at "
1267  "least two-dimensional");
1268  dir = 1;
1269  break;
1271  ASSERTL1(m_geom->GetCoordim() == 3,
1272  "Cannot call eWeakDeriv2 in a "
1273  "coordinate system which is not "
1274  "three-dimensional");
1275  dir = 2;
1276  break;
1277  default:
1278  break;
1279  }
1280 
1281  MatrixKey deriv0key(StdRegions::eWeakDeriv0,
1282  mkey.GetShapeType(), *this);
1283 
1284  DNekMatSharedPtr WeakDerivStd = GetStdMatrix(deriv0key);
1285  fac = m_metricinfo->GetDerivFactors(ptsKeys)[dir][0]*
1286  m_metricinfo->GetJac(ptsKeys)[0];
1287 
1288  returnval = MemoryManager<DNekScalMat>::
1289  AllocateSharedPtr(fac,WeakDerivStd);
1290  }
1291  }
1292  break;
1294  {
1295  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1296  {
1297  fac = 1.0;
1298  goto UseLocRegionsMatrix;
1299  }
1300  else
1301  {
1302  int coordim = m_geom->GetCoordim();
1303  fac = 0.0;
1304  for (int i = 0; i < coordim; ++i)
1305  {
1306  fac += m_metricinfo->GetDerivFactors(ptsKeys)[i][0]*
1307  m_metricinfo->GetDerivFactors(ptsKeys)[i][0];
1308  }
1309  fac *= m_metricinfo->GetJac(ptsKeys)[0];
1310  goto UseStdRegionsMatrix;
1311  }
1312  }
1313  break;
1315  {
1316  NekDouble factor =
1317  mkey.GetConstFactor(StdRegions::eFactorLambda);
1318  MatrixKey masskey(StdRegions::eMass,
1319  mkey.GetShapeType(), *this);
1320  DNekScalMat &MassMat = *(this->m_matrixManager[masskey]);
1321  MatrixKey lapkey(StdRegions::eLaplacian, mkey.GetShapeType(),
1322  *this, mkey.GetConstFactors(),
1323  mkey.GetVarCoeffs());
1324  DNekScalMat &LapMat = *(this->m_matrixManager[lapkey]);
1325 
1326  int rows = LapMat.GetRows();
1327  int cols = LapMat.GetColumns();
1328 
1329  DNekMatSharedPtr helm =
1331 
1332  NekDouble one = 1.0;
1333  (*helm) = LapMat + factor*MassMat;
1334 
1335  returnval =
1337  }
1338  break;
1343  {
1344  NekDouble one = 1.0;
1345 
1346  DNekMatSharedPtr mat = GenMatrix(mkey);
1347  returnval =
1349  }
1350  break;
1352  {
1353  NekDouble one = 1.0;
1354 
1355 // StdRegions::StdMatrixKey hkey(StdRegions::eHybridDGHelmholtz,
1356 // DetShapeType(),*this,
1357 // mkey.GetConstant(0),
1358 // mkey.GetConstant(1));
1359  MatrixKey hkey(StdRegions::eHybridDGHelmholtz,
1360  DetShapeType(),
1361  *this, mkey.GetConstFactors(),
1362  mkey.GetVarCoeffs());
1363  DNekMatSharedPtr mat = GenMatrix(hkey);
1364 
1365  mat->Invert();
1366  returnval =
1368  }
1369  break;
1371  {
1372  DNekMatSharedPtr m_Ix;
1373  Array<OneD, NekDouble> coords(1, 0.0);
1374  StdRegions::ConstFactorMap factors = mkey.GetConstFactors();
1375  int vertex = (int)factors[StdRegions::eFactorGaussVertex];
1376 
1377  coords[0] = (vertex == 0) ? -1.0 : 1.0;
1378 
1379  m_Ix = m_base[0]->GetI(coords);
1380  returnval =
1382  }
1383  break;
1384 
1385  UseLocRegionsMatrix:
1386  {
1387  DNekMatSharedPtr mat = GenMatrix(mkey);
1388  returnval =
1390  }
1391  break;
1392  UseStdRegionsMatrix:
1393  {
1394  DNekMatSharedPtr mat = GetStdMatrix(mkey);
1395  returnval =
1397  }
1398  break;
1399  default:
1400  NEKERROR(ErrorUtil::efatal, "Matrix creation not defined");
1401  break;
1402  }
1403 
1404  return returnval;
1405  }
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:158
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:220
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:126
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:251
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:125
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:700
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
double NekDouble
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:213
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
Array< OneD, LibUtilities::BasisSharedPtr > m_base
Geometry is curved or has non-constant factors.
NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag > DNekScalMat
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: SegExp.h:246
DNekScalBlkMatSharedPtr Nektar::LocalRegions::SegExp::CreateStaticCondMatrix ( const MatrixKey mkey)
protected
Todo:
Really need a constructor which takes Arrays

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

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

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

1623  {
1624  // get Mass matrix inverse
1625  MatrixKey masskey(StdRegions::eInvMass,
1626  DetShapeType(),*this);
1627  DNekScalMatSharedPtr matsys = m_matrixManager[masskey];
1628 
1629  NekVector<NekDouble> in(m_ncoeffs,inarray,eCopy);
1630  NekVector<NekDouble> out(m_ncoeffs,outarray,eWrapper);
1631 
1632  out = (*matsys)*in;
1633  }
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:246
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 1578 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().

1581  {
1582 
1583  int m;
1584  NekDouble sgn = 1;
1585 
1586  ASSERTL1(&inarray[0] != &outarray[0],
1587  "inarray and outarray can not be the same");
1588  switch(GetBasisType(0))
1589  {
1591  //Swap vertices
1592  outarray[0] = inarray[1];
1593  outarray[1] = inarray[0];
1594  // negate odd modes
1595  for(m = 2; m < m_ncoeffs; ++m)
1596  {
1597  outarray[m] = sgn*inarray[m];
1598  sgn = -sgn;
1599  }
1600  break;
1603  for(m = 0; m < m_ncoeffs; ++m)
1604  {
1605  outarray[m_ncoeffs-1-m] = inarray[m];
1606  }
1607  break;
1608  default:
1609  ASSERTL0(false,"This basis is not allowed in this method");
1610  break;
1611  }
1612  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
Principle Modified Functions .
Definition: BasisType.h:49
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:54
double NekDouble
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
Lagrange for SEM basis .
Definition: BasisType.h:53
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
void Nektar::LocalRegions::SegExp::v_ComputeVertexNormal ( const int  vertex)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 884 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, Vmath::Smul(), and Nektar::NekMeshUtils::vert.

885  {
886  int i;
887  const SpatialDomains::GeomFactorsSharedPtr &geomFactors =
888  GetGeom()->GetMetricInfo();
889  SpatialDomains::GeomType type = geomFactors->GetGtype();
890  const Array<TwoD, const NekDouble> &gmat =
891  geomFactors->GetDerivFactors(GetPointsKeys());
892  int nqe = 1;
893  int vCoordDim = GetCoordim();
894 
895  m_vertexNormals[vertex] =
896  Array<OneD, Array<OneD, NekDouble> >(vCoordDim);
897  Array<OneD, Array<OneD, NekDouble> > &normal =
898  m_vertexNormals[vertex];
899  for (i = 0; i < vCoordDim; ++i)
900  {
901  normal[i] = Array<OneD, NekDouble>(nqe);
902  }
903 
904  // Regular geometry case
905  if ((type == SpatialDomains::eRegular) ||
907  {
908  NekDouble vert;
909  // Set up normals
910  switch (vertex)
911  {
912  case 0:
913  for(i = 0; i < vCoordDim; ++i)
914  {
915  Vmath::Fill(nqe, -gmat[i][0], normal[i], 1);
916  }
917  break;
918  case 1:
919  for(i = 0; i < vCoordDim; ++i)
920  {
921  Vmath::Fill(nqe, gmat[i][0], normal[i], 1);
922  }
923  break;
924  default:
925  ASSERTL0(false,
926  "point is out of range (point < 2)");
927  }
928 
929  // normalise
930  vert = 0.0;
931  for (i =0 ; i < vCoordDim; ++i)
932  {
933  vert += normal[i][0]*normal[i][0];
934  }
935  vert = 1.0/sqrt(vert);
936  for (i = 0; i < vCoordDim; ++i)
937  {
938  Vmath::Smul(nqe, vert, normal[i], 1, normal[i], 1);
939  }
940  }
941  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
std::map< int, NormalVector > m_vertexNormals
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:46
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
double NekDouble
SpatialDomains::GeometrySharedPtr GetGeom() const
Definition: Expansion.cpp:148
boost::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition: GeomFactors.h:62
Geometry is straight-sided with constant geometric factors.
GeomType
Indicates the type of element geometry.
DNekMatSharedPtr Nektar::LocalRegions::SegExp::v_CreateStdMatrix ( const StdRegions::StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 1186 of file SegExp.cpp.

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

1188  {
1189  LibUtilities::BasisKey bkey = m_base[0]->GetBasisKey();
1192 
1193  return tmp->GetStdMatrix(mkey);
1194  }
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 1175 of file SegExp.cpp.

References m_staticCondMatrixManager.

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

Unpack data from input file assuming it comes from.

Reimplemented from Nektar::StdRegions::StdExpansion.

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

845  {
846  switch(m_base[0]->GetBasisType())
847  {
849  {
850  int fillorder = min((int) nummodes[mode_offset],m_ncoeffs);
851 
852  Vmath::Zero(m_ncoeffs,coeffs,1);
853  Vmath::Vcopy(fillorder,&data[0],1,&coeffs[0],1);
854  }
855  break;
857  {
858  // Assume that input is also Gll_Lagrange
859  // but no way to check;
860  LibUtilities::PointsKey p0(
861  nummodes[mode_offset],
864  p0,data, m_base[0]->GetPointsKey(), coeffs);
865  }
866  break;
868  {
869  // Assume that input is also Gauss_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;
878  default:
879  ASSERTL0(false,
880  "basis is either not set up or not hierarchicial");
881  }
882  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
Principle Modified Functions .
Definition: BasisType.h:49
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:54
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:47
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
void Interp1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
this function interpolates a 1D function evaluated at the quadrature points of the basis fbasis0 to ...
Definition: Interp.cpp:54
Lagrange for SEM basis .
Definition: BasisType.h:53
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:359
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:50
void Nektar::LocalRegions::SegExp::v_FwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

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

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

Inputs:

  • inarray: array of physical quadrature points to be transformed

Outputs:

  • outarray: updated array of expansion coefficients.

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 375 of file SegExp.cpp.

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

Referenced by v_FwdTrans_BndConstrained().

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

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 399 of file SegExp.cpp.

References ASSERTL0, 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().

402  {
403  if(m_base[0]->Collocation())
404  {
405  Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
406  }
407  else
408  {
409  int nInteriorDofs = m_ncoeffs-2;
410  int offset;
411 
412  switch (m_base[0]->GetBasisType())
413  {
415  {
416  offset = 1;
417  }
418  break;
420  {
421  nInteriorDofs = m_ncoeffs;
422  offset = 0;
423  }
424  break;
427  {
428  ASSERTL1(m_base[0]->GetPointsType() == LibUtilities::eGaussLobattoLegendre,"Cannot use FwdTrans_BndConstrained method with non GLL points");
429  offset = 2;
430  }
431  break;
432  default:
433  ASSERTL0(false,"This type of FwdTrans is not defined"
434  "for this expansion type");
435  }
436 
437  fill(outarray.get(), outarray.get()+m_ncoeffs, 0.0 );
438 
440  {
441 
442  outarray[GetVertexMap(0)] = inarray[0];
443  outarray[GetVertexMap(1)] =
444  inarray[m_base[0]->GetNumPoints()-1];
445 
446  if (m_ncoeffs>2)
447  {
448  // ideally, we would like to have tmp0 to be replaced
449  // by outarray (currently MassMatrixOp does not allow
450  // aliasing)
451  Array<OneD, NekDouble> tmp0(m_ncoeffs);
452  Array<OneD, NekDouble> tmp1(m_ncoeffs);
453 
454  StdRegions::StdMatrixKey stdmasskey(
456  MassMatrixOp(outarray,tmp0,stdmasskey);
457  v_IProductWRTBase(inarray,tmp1);
458 
459  Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
460 
461  // get Mass matrix inverse (only of interior DOF)
462  MatrixKey masskey(
464  DNekScalMatSharedPtr matsys =
465  (m_staticCondMatrixManager[masskey])->GetBlock(1,1);
466 
467  Blas::Dgemv('N',nInteriorDofs,nInteriorDofs,
468  matsys->Scale(),
469  &((matsys->GetOwnedMatrix())->GetPtr())[0],
470  nInteriorDofs,tmp1.get()+offset,1,0.0,
471  outarray.get()+offset,1);
472  }
473  }
474  else
475  {
476  SegExp::v_FwdTrans(inarray, outarray);
477 
478  }
479  }
480  }
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:470
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:971
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: SegExp.h:248
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Forward transform from physical quadrature space stored in inarray and evaluate the expansion coeffic...
Definition: SegExp.cpp:375
Principle Modified Functions .
Definition: BasisType.h:49
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:54
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Inner product of inarray over region with respect to the expansion basis (this)->_Base[0] and return ...
Definition: SegExp.cpp:502
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
int GetVertexMap(const int localVertexId, bool useCoeffPacking=false)
Definition: StdExpansion.h:826
Principle Modified Functions .
Definition: BasisType.h:50
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:165
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:329
Lagrange for SEM basis .
Definition: BasisType.h:53
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:191
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:50
DNekMatSharedPtr Nektar::LocalRegions::SegExp::v_GenMatrix ( const StdRegions::StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

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

1409  {
1410  DNekMatSharedPtr returnval;
1411 
1412  switch (mkey.GetMatrixType())
1413  {
1420  returnval = Expansion1D::v_GenMatrix(mkey);
1421  break;
1422  default:
1423  returnval = StdSegExp::v_GenMatrix(mkey);
1424  break;
1425  }
1426 
1427  return returnval;
1428  }
virtual DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey)
Definition: Expansion1D.cpp:43
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
const LibUtilities::BasisSharedPtr & Nektar::LocalRegions::SegExp::v_GetBasis ( int  dir) const
protectedvirtual

Definition at line 820 of file SegExp.cpp.

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

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

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

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

Reimplemented from Nektar::StdRegions::StdExpansion1D.

Definition at line 794 of file SegExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 694 of file SegExp.cpp.

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

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

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1180 of file SegExp.cpp.

References m_matrixManager.

1181  {
1182  return m_matrixManager[mkey];
1183  }
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: SegExp.h:246
DNekScalBlkMatSharedPtr Nektar::LocalRegions::SegExp::v_GetLocStaticCondMatrix ( const MatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1169 of file SegExp.cpp.

References m_staticCondMatrixManager.

1171  {
1172  return m_staticCondMatrixManager[mkey];
1173  }
LibUtilities::NekManager< MatrixKey, DNekScalBlkMat, MatrixKey::opLess > m_staticCondMatrixManager
Definition: SegExp.h:248
int Nektar::LocalRegions::SegExp::v_GetNcoeffs ( void  ) const
protectedvirtual

Definition at line 815 of file SegExp.cpp.

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

816  {
817  return m_ncoeffs;
818  }
int Nektar::LocalRegions::SegExp::v_GetNumPoints ( const int  dir) const
protectedvirtual

Definition at line 810 of file SegExp.cpp.

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

811  {
812  return GetNumPoints(dir);
813  }
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 799 of file SegExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 782 of file SegExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 788 of file SegExp.cpp.

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

789  {
791  ::AllocateSharedPtr(m_base[0]->GetBasisKey());
792  }
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 739 of file SegExp.cpp.

References v_GetVertexPhysVals().

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

707  {
708  int nquad = m_base[0]->GetNumPoints();
709 
711  {
712  switch (vertex)
713  {
714  case 0:
715  outarray = inarray[0];
716  break;
717  case 1:
718  outarray = inarray[nquad - 1];
719  break;
720  }
721  }
722  else
723  {
725  factors[StdRegions::eFactorGaussVertex] = vertex;
726 
727  StdRegions::StdMatrixKey key(
729  DetShapeType(),*this,factors);
730 
731  DNekScalMatSharedPtr mat_gauss = m_matrixManager[key];
732 
733  outarray = Blas::Ddot(nquad, mat_gauss->GetOwnedMatrix()
734  ->GetPtr().get(), 1, &inarray[0], 1);
735  }
736  }
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:251
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:47
T Ddot(int n, const Array< OneD, const T > &w, const int incw, const Array< OneD, const T > &x, const int incx, const Array< OneD, const int > &y, const int incy)
Definition: VmathArray.hpp:434
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:216
Array< OneD, LibUtilities::BasisSharedPtr > m_base
LibUtilities::NekManager< MatrixKey, DNekScalMat, MatrixKey::opLess > m_matrixManager
Definition: SegExp.h:246
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 1054 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().

1058  {
1059  int nquad = m_base[0]->GetNumPoints();
1060  const Array<TwoD, const NekDouble>& gmat =
1061  m_metricinfo->GetDerivFactors(GetPointsKeys());
1062  const NekDouble lambda =
1063  mkey.GetConstFactor(StdRegions::eFactorLambda);
1064 
1065  Array<OneD,NekDouble> physValues(nquad);
1066  Array<OneD,NekDouble> dPhysValuesdx(nquad);
1067  Array<OneD,NekDouble> wsp(m_ncoeffs);
1068 
1069  BwdTrans(inarray, physValues);
1070 
1071  // mass matrix operation
1072  v_IProductWRTBase((m_base[0]->GetBdata()),physValues,wsp,1);
1073 
1074  // Laplacian matrix operation
1075  switch (m_geom->GetCoordim())
1076  {
1077  case 1:
1078  {
1079  PhysDeriv(physValues,dPhysValuesdx);
1080 
1081  // multiply with the proper geometric factors
1082  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1083  {
1084  Vmath::Vmul(nquad,
1085  &gmat[0][0],1,dPhysValuesdx.get(),1,
1086  dPhysValuesdx.get(),1);
1087  }
1088  else
1089  {
1090  Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.get(), 1);
1091  }
1092  }
1093  break;
1094  case 2:
1095  {
1096  Array<OneD,NekDouble> dPhysValuesdy(nquad);
1097 
1098  PhysDeriv(physValues, dPhysValuesdx, dPhysValuesdy);
1099 
1100  // multiply with the proper geometric factors
1101  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1102  {
1103  Vmath::Vmul (nquad,
1104  &gmat[0][0], 1, dPhysValuesdx.get(), 1,
1105  dPhysValuesdx.get(), 1);
1106  Vmath::Vvtvp(nquad,
1107  &gmat[1][0], 1, dPhysValuesdy.get(), 1,
1108  dPhysValuesdx.get(), 1,
1109  dPhysValuesdx.get(), 1);
1110  }
1111  else
1112  {
1113  Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.get(), 1);
1114  Blas::Daxpy(nquad,
1115  gmat[1][0], dPhysValuesdy.get(), 1,
1116  dPhysValuesdx.get(), 1);
1117  }
1118  }
1119  break;
1120  case 3:
1121  {
1122  Array<OneD,NekDouble> dPhysValuesdy(nquad);
1123  Array<OneD,NekDouble> dPhysValuesdz(nquad);
1124 
1125  PhysDeriv(physValues, dPhysValuesdx,
1126  dPhysValuesdy, dPhysValuesdz);
1127 
1128  // multiply with the proper geometric factors
1129  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1130  {
1131  Vmath::Vmul (nquad,
1132  &gmat[0][0], 1, dPhysValuesdx.get(), 1,
1133  dPhysValuesdx.get(), 1);
1134  Vmath::Vvtvp(nquad,
1135  &gmat[1][0], 1, dPhysValuesdy.get(), 1,
1136  dPhysValuesdx.get(), 1,
1137  dPhysValuesdx.get(), 1);
1138  Vmath::Vvtvp(nquad,
1139  &gmat[2][0], 1, dPhysValuesdz.get(), 1,
1140  dPhysValuesdx.get(), 1,
1141  dPhysValuesdx.get(), 1);
1142  }
1143  else
1144  {
1145  Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.get(), 1);
1146  Blas::Daxpy(nquad,
1147  gmat[1][0], dPhysValuesdy.get(), 1,
1148  dPhysValuesdx.get(), 1);
1149  Blas::Daxpy(nquad,
1150  gmat[2][0], dPhysValuesdz.get(),
1151  1, dPhysValuesdx.get(), 1);
1152  }
1153  }
1154  break;
1155  default:
1156  ASSERTL0(false,"Wrong number of dimensions");
1157  break;
1158  }
1159 
1160  v_IProductWRTBase(m_base[0]->GetDbdata(),dPhysValuesdx,outarray,1);
1161  Blas::Daxpy(m_ncoeffs, lambda, wsp.get(), 1, outarray.get(), 1);
1162  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:428
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:126
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:125
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Inner product of inarray over region with respect to the expansion basis (this)->_Base[0] and return ...
Definition: SegExp.cpp:502
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:525
Array< OneD, LibUtilities::BasisSharedPtr > m_base
Geometry is curved or has non-constant factors.
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:169
NekDouble Nektar::LocalRegions::SegExp::v_Integral ( const Array< OneD, const NekDouble > &  inarray)
protectedvirtual

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

Inputs:

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

Outputs:

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

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 116 of file SegExp.cpp.

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

118  {
119  int nquad0 = m_base[0]->GetNumPoints();
120  Array<OneD, const NekDouble> jac = m_metricinfo->GetJac(GetPointsKeys());
121  NekDouble ival;
122  Array<OneD,NekDouble> tmp(nquad0);
123 
124  // multiply inarray with Jacobian
125  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
126  {
127  Vmath::Vmul(nquad0, jac, 1, inarray, 1, tmp,1);
128  }
129  else
130  {
131  Vmath::Smul(nquad0, jac[0], inarray, 1, tmp, 1);
132  }
133 
134  // call StdSegExp version;
135  ival = StdSegExp::v_Integral(tmp);
136  //ival = StdSegExp::Integral(tmp);
137  return ival;
138  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:126
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
double NekDouble
Array< OneD, LibUtilities::BasisSharedPtr > m_base
Geometry is curved or has non-constant factors.
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:169
void Nektar::LocalRegions::SegExp::v_IProductWRTBase ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

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

Wrapper call to SegExp::IProduct_WRT_B

Input:

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

Output:

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

Reimplemented from Nektar::StdRegions::StdSegExp.

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

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

542  {
543  int nquad0 = m_base[0]->GetNumPoints();
544  Array<OneD, const NekDouble> jac = m_metricinfo->GetJac(GetPointsKeys());
545  Array<OneD,NekDouble> tmp(nquad0);
546 
547 
548  // multiply inarray with Jacobian
549  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
550  {
551  Vmath::Vmul(nquad0, jac, 1, inarray, 1, tmp, 1);
552  }
553  else
554  {
555  Vmath::Smul(nquad0, jac[0], inarray, 1, tmp, 1);
556  }
557  StdSegExp::v_IProductWRTBase(base,tmp,outarray,coll_check);
558  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:126
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
Array< OneD, LibUtilities::BasisSharedPtr > m_base
Geometry is curved or has non-constant factors.
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:169
void Nektar::LocalRegions::SegExp::v_IProductWRTDerivBase ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

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

565  {
566  int nquad = m_base[0]->GetNumPoints();
567  const Array<TwoD, const NekDouble>& gmat =
568  m_metricinfo->GetDerivFactors(GetPointsKeys());
569 
570  Array<OneD, NekDouble> tmp1(nquad);
571 
572  switch(dir)
573  {
574  case 0:
575  {
576  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
577  {
578  Vmath::Vmul(nquad,gmat[0],1,inarray,1,tmp1,1);
579  }
580  else
581  {
582  Vmath::Smul(nquad, gmat[0][0], inarray, 1, tmp1, 1);
583  }
584  }
585  break;
586  case 1:
587  {
588  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
589  {
590  Vmath::Vmul(nquad,gmat[1],1,inarray,1,tmp1,1);
591  }
592  else
593  {
594  Vmath::Smul(nquad, gmat[1][0], inarray, 1, tmp1, 1);
595  }
596  }
597  break;
598  case 2:
599  {
600  ASSERTL1(m_geom->GetCoordim() == 3,"input dir is out of range");
601  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
602  {
603  Vmath::Vmul(nquad,gmat[2],1,inarray,1,tmp1,1);
604  }
605  else
606  {
607  Vmath::Smul(nquad, gmat[2][0], inarray, 1, tmp1, 1);
608  }
609  }
610  break;
611  default:
612  {
613  ASSERTL1(false,"input dir is out of range");
614  }
615  break;
616  }
617  v_IProductWRTBase(m_base[0]->GetDbdata(),tmp1,outarray,1);
618  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:126
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:125
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Inner product of inarray over region with respect to the expansion basis (this)->_Base[0] and return ...
Definition: SegExp.cpp:502
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
Array< OneD, LibUtilities::BasisSharedPtr > m_base
Geometry is curved or has non-constant factors.
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:169
void Nektar::LocalRegions::SegExp::v_LaplacianMatrixOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdRegions::StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

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

952  {
953  int nquad = m_base[0]->GetNumPoints();
954  const Array<TwoD, const NekDouble>& gmat =
955  m_metricinfo->GetDerivFactors(GetPointsKeys());
956 
957  Array<OneD,NekDouble> physValues(nquad);
958  Array<OneD,NekDouble> dPhysValuesdx(nquad);
959 
960  BwdTrans(inarray,physValues);
961 
962  // Laplacian matrix operation
963  switch (m_geom->GetCoordim())
964  {
965  case 1:
966  {
967  PhysDeriv(physValues,dPhysValuesdx);
968 
969  // multiply with the proper geometric factors
970  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
971  {
972  Vmath::Vmul(nquad,
973  &gmat[0][0],1,dPhysValuesdx.get(),1,
974  dPhysValuesdx.get(),1);
975  }
976  else
977  {
978  Blas::Dscal(nquad,
979  gmat[0][0], dPhysValuesdx.get(), 1);
980  }
981  }
982  break;
983  case 2:
984  {
985  Array<OneD,NekDouble> dPhysValuesdy(nquad);
986 
987  PhysDeriv(physValues,dPhysValuesdx,dPhysValuesdy);
988 
989  // multiply with the proper geometric factors
990  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
991  {
992  Vmath::Vmul (nquad,
993  &gmat[0][0],1,dPhysValuesdx.get(),1,
994  dPhysValuesdx.get(),1);
995  Vmath::Vvtvp(nquad,
996  &gmat[1][0],1,dPhysValuesdy.get(),1,
997  dPhysValuesdx.get(),1,
998  dPhysValuesdx.get(),1);
999  }
1000  else
1001  {
1002  Blas::Dscal(nquad,
1003  gmat[0][0], dPhysValuesdx.get(), 1);
1004  Blas::Daxpy(nquad,
1005  gmat[1][0], dPhysValuesdy.get(), 1,
1006  dPhysValuesdx.get(), 1);
1007  }
1008  }
1009  break;
1010  case 3:
1011  {
1012  Array<OneD,NekDouble> dPhysValuesdy(nquad);
1013  Array<OneD,NekDouble> dPhysValuesdz(nquad);
1014 
1015  PhysDeriv(physValues,dPhysValuesdx,
1016  dPhysValuesdy,dPhysValuesdz);
1017 
1018  // multiply with the proper geometric factors
1019  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1020  {
1021  Vmath::Vmul (nquad,
1022  &gmat[0][0], 1, dPhysValuesdx.get(), 1,
1023  dPhysValuesdx.get(),1);
1024  Vmath::Vvtvp(nquad,
1025  &gmat[1][0], 1, dPhysValuesdy.get(), 1,
1026  dPhysValuesdx.get(),1,
1027  dPhysValuesdx.get(),1);
1028  Vmath::Vvtvp(nquad,
1029  &gmat[2][0], 1, dPhysValuesdz.get(), 1,
1030  dPhysValuesdx.get(),1,
1031  dPhysValuesdx.get(),1);
1032  }
1033  else
1034  {
1035  Blas::Dscal(nquad, gmat[0][0], dPhysValuesdx.get(), 1);
1036  Blas::Daxpy(nquad,
1037  gmat[1][0], dPhysValuesdy.get(), 1,
1038  dPhysValuesdx.get(), 1);
1039  Blas::Daxpy(nquad,
1040  gmat[2][0], dPhysValuesdz.get(), 1,
1041  dPhysValuesdx.get(), 1);
1042  }
1043  }
1044  break;
1045  default:
1046  ASSERTL0(false,"Wrong number of dimensions");
1047  break;
1048  }
1049 
1050  v_IProductWRTBase(m_base[0]->GetDbdata(),dPhysValuesdx,outarray,1);
1051  }
const LibUtilities::PointsKeyVector GetPointsKeys() const
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:428
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:126
SpatialDomains::GeometrySharedPtr m_geom
Definition: Expansion.h:125
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Inner product of inarray over region with respect to the expansion basis (this)->_Base[0] and return ...
Definition: SegExp.cpp:502
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:525
Array< OneD, LibUtilities::BasisSharedPtr > m_base
Geometry is curved or has non-constant factors.
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:169
SpatialDomains::GeomType Nektar::LocalRegions::SegExp::v_MetricInfoType ( )
protectedvirtual

Definition at line 805 of file SegExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

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

624  {
625  int nq = m_base[0]->GetNumPoints();
626  Array<OneD, NekDouble > Fn(nq);
627 // cout << "I am segment " << GetGeom()->GetGlobalID() << endl;
628 // cout << "I want edge " << GetLeftAdjacentElementEdge() << endl;
629 // @TODO: This routine no longer makes sense as a normal is not unique to an edge
630  const Array<OneD, const Array<OneD, NekDouble> >
631  &normals =
634  Vmath::Vmul (nq, &Fx[0], 1, &normals[0][0], 1, &Fn[0], 1);
635  Vmath::Vvtvp(nq, &Fy[0], 1, &normals[1][0], 1, &Fn[0], 1, &Fn[0], 1);
636 
637  v_IProductWRTBase(Fn,outarray);
638  }
const NormalVector & GetEdgeNormal(const int edge) const
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:428
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Inner product of inarray over region with respect to the expansion basis (this)->_Base[0] and return ...
Definition: SegExp.cpp:502
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:169
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 640 of file SegExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 826 of file SegExp.cpp.

827  {
828  return 2;
829  }
int Nektar::LocalRegions::SegExp::v_NumDGBndryCoeffs ( ) const
protectedvirtual

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 832 of file SegExp.cpp.

833  {
834  return 2;
835  }
void Nektar::LocalRegions::SegExp::v_PhysDeriv ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d0,
Array< OneD, NekDouble > &  out_d1 = NullNekDouble1DArray,
Array< OneD, NekDouble > &  out_d2 = NullNekDouble1DArray 
)
protectedvirtual

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

This is a wrapper around StdExpansion1D::Tensor_Deriv

Input:

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

Output:

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

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 164 of file SegExp.cpp.

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

Referenced by v_PhysDeriv_n().

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

Calculate the derivative of the physical points in a given direction.

See also
StdRegions::StdExpansion::PhysDeriv

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 317 of file SegExp.cpp.

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

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

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

Parameters
inarrayfunction to derive
out_dnresult of the derivative operation

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 267 of file SegExp.cpp.

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

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

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

Parameters
inarrayfunction to derive
out_dsresult of the derivative operation

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 227 of file SegExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdSegExp.

Definition at line 665 of file SegExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 755 of file SegExp.cpp.

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

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

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

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

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