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::QuadExp Class Reference

#include <QuadExp.h>

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

Public Member Functions

 QuadExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const SpatialDomains::QuadGeomSharedPtr &geom)
 Constructor using BasisKey class for quadrature points and order definition.
 QuadExp (const QuadExp &T)
virtual ~QuadExp ()
- Public Member Functions inherited from Nektar::StdRegions::StdQuadExp
 StdQuadExp ()
 StdQuadExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb)
 Constructor using BasisKey class for quadrature points and order definition.
 StdQuadExp (const StdQuadExp &T)
 Copy Constructor.
 ~StdQuadExp ()
 Destructor.
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion2D
 StdExpansion2D ()
 StdExpansion2D (int numcoeffs, const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb)
 StdExpansion2D (const StdExpansion2D &T)
virtual ~StdExpansion2D ()
void PhysTensorDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d0, Array< OneD, NekDouble > &outarray_d1)
 Calculate the 2D derivative in the local tensor/collapsed coordinate at the physical points.
NekDouble Integral (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &w0, const Array< OneD, const NekDouble > &w1)
void BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
void IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion
 StdExpansion ()
 Default Constructor.
 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.
 StdExpansion (const StdExpansion &T)
 Copy Constructor.
virtual ~StdExpansion ()
 Destructor.
int GetNumBases () const
 This function returns the number of 1D bases used in the expansion.
const Array< OneD, const
LibUtilities::BasisSharedPtr > & 
GetBase () const
 This function gets the shared point to basis.
const
LibUtilities::BasisSharedPtr
GetBasis (int dir) const
 This function gets the shared point to basis in the dir direction.
int GetNcoeffs (void) const
 This function returns the total number of coefficients used in the expansion.
int GetTotPoints () const
 This function returns the total number of quadrature points used in the element.
LibUtilities::BasisType GetBasisType (const int dir) const
 This function returns the type of basis used in the dir direction.
int GetBasisNumModes (const int dir) const
 This function returns the number of expansion modes in the dir direction.
int EvalBasisNumModesMax (void) const
 This function returns the maximum number of expansion modes over all local directions.
LibUtilities::PointsType GetPointsType (const int dir) const
 This function returns the type of quadrature points used in the dir direction.
int GetNumPoints (const int dir) const
 This function returns the number of quadrature points in the dir direction.
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.
int GetNverts () const
 This function returns the number of vertices of the expansion domain.
int GetNedges () const
 This function returns the number of edges of the expansion domain.
int GetEdgeNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th edge.
int GetTotalEdgeIntNcoeffs () const
int GetEdgeNumPoints (const int i) const
 This function returns the number of quadrature points belonging to the i-th edge.
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.
int GetFaceNcoeffs (const int i) const
 This function returns the number of expansion coefficients belonging to the i-th face.
int GetFaceIntNcoeffs (const int i) const
int GetTotalFaceIntNcoeffs () const
LibUtilities::PointsKey GetFacePointsKey (const int i, const int j) const
int NumBndryCoeffs (void) const
int NumDGBndryCoeffs (void) const
LibUtilities::BasisType GetEdgeBasisType (const int i) const
 This function returns the type of expansion basis on the i-th edge.
int GetNfaces () const
 This function returns the number of faces of the expansion domain.
int GetNtrace () const
 Returns the number of trace elements connected to this element.
LibUtilities::ShapeType DetShapeType () const
 This function returns the shape of the expansion domain.
int GetShapeDimension () const
bool IsBoundaryInteriorExpansion ()
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Backward transformation from coefficient space to physical space.
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs the Forward transformation from physical space to coefficient space.
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.
void FillMode (const int mode, Array< OneD, NekDouble > &outarray)
 This function fills the array outarray with the mode-th mode of the expansion.
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
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.
void SetElmtId (const int id)
 Set the element id of this expansion when used in a list by returning value of m_elmt_id.
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
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
DNekMatSharedPtr GetStdMatrix (const StdMatrixKey &mkey)
DNekBlkMatSharedPtr GetStdStaticCondMatrix (const StdMatrixKey &mkey)
IndexMapValuesSharedPtr GetIndexMap (const IndexMapKey &ikey)
const Array< OneD, const
NekDouble > & 
GetPhysNormals (void)
void SetPhysNormals (Array< OneD, const NekDouble > &normal)
virtual void SetUpPhysNormals (const int edge)
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
void NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
DNekScalBlkMatSharedPtr GetLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
void DropLocStaticCondMatrix (const LocalRegions::MatrixKey &mkey)
StdRegions::Orientation GetForient (int face)
StdRegions::Orientation GetEorient (int edge)
StdRegions::Orientation GetPorient (int point)
StdRegions::Orientation GetCartesianEorient (int edge)
void SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
void SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
int CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
void ExtractDataToCoeffs (const NekDouble *data, const std::vector< unsigned int > &nummodes, const int nmodes_offset, NekDouble *coeffs)
NekDouble StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
int GetCoordim ()
void GetBoundaryMap (Array< OneD, unsigned int > &outarray)
void GetInteriorMap (Array< OneD, unsigned int > &outarray)
int GetVertexMap (const int localVertexId, bool useCoeffPacking=false)
void GetEdgeInteriorMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
void GetFaceInteriorMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
void GetEdgeToElementMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
void GetFaceToElementMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int nummodesA=-1, int nummodesB=-1)
void GetEdgePhysVals (const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Extract the physical values along edge edge from inarray into outarray following the local edge orientation and point distribution defined by defined in EdgeExp.
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).
void GetFacePhysVals (const int face, const boost::shared_ptr< StdExpansion > &FaceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=eNoOrientation)
void MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
DNekMatSharedPtr CreateGeneralMatrix (const StdMatrixKey &mkey)
 this function generates the mass matrix $\mathbf{M}[i][j] = \int \phi_i(\mathbf{x}) \phi_j(\mathbf{x}) d\mathbf{x}$
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.
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.
void LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
 Convert local cartesian coordinate xi into local collapsed coordinates eta.
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.
virtual const Array< OneD,
const NekDouble > & 
v_GetPhysNormals (void)
virtual void v_SetPhysNormals (Array< OneD, const NekDouble > &normal)
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
virtual StdRegions::Orientation v_GetForient (int face)
virtual StdRegions::Orientation v_GetPorient (int point)
NekDouble Linf (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete $ L_\infty$ error $ |\epsilon|_\infty = \max |u - u_{exact}|$ where $ u_{exact}$ is given by the array sol.
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.
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.
const NormalVectorGetEdgeNormal (const int edge) const
void ComputeEdgeNormal (const int edge)
void NegateEdgeNormal (const int edge)
bool EdgeNormalNegated (const int edge)
void ComputeFaceNormal (const int face)
void NegateFaceNormal (const int face)
void ComputeVertexNormal (const int vertex)
const NormalVectorGetFaceNormal (const int face) const
const NormalVectorGetVertexNormal (const int vertex) const
const NormalVectorGetSurfaceNormal (const int id) const
const LibUtilities::PointsKeyVector GetPointsKeys () const
Array< OneD, unsigned int > GetEdgeInverseBoundaryMap (int eid)
Array< OneD, unsigned int > GetFaceInverseBoundaryMap (int fid, StdRegions::Orientation faceOrient=eNoOrientation)
DNekMatSharedPtr BuildInverseTransformationMatrix (const DNekScalMatSharedPtr &m_transformationmatrix)
void PhysInterpToSimplexEquiSpaced (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function performs an interpolation from the physical space points provided at input into an array of equispaced points which are not the collapsed coordinate. So for a tetrahedron you will only get a tetrahedral number of values.
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.
template<class T >
boost::shared_ptr< T > as ()
- Public Member Functions inherited from Nektar::LocalRegions::Expansion2D
 Expansion2D (SpatialDomains::Geometry2DSharedPtr pGeom)
virtual ~Expansion2D ()
void SetTraceToGeomOrientation (Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, NekDouble > &inout)
ExpansionSharedPtr GetEdgeExp (int edge, bool SetUpNormal=true)
void SetEdgeExp (const int edge, ExpansionSharedPtr &e)
void AddNormTraceInt (const int dir, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &edgeCoeffs, Array< OneD, NekDouble > &outarray)
void AddNormTraceInt (const int dir, Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, NekDouble > &outarray, const StdRegions::VarCoeffMap &varcoeffs)
void AddEdgeBoundaryInt (const int edge, ExpansionSharedPtr &EdgeExp, Array< OneD, NekDouble > &edgePhys, Array< OneD, NekDouble > &outarray, const StdRegions::VarCoeffMap &varcoeffs=StdRegions::NullVarCoeffMap)
void AddHDGHelmholtzEdgeTerms (const NekDouble tau, const int edge, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, NekDouble > &edgePhys, const StdRegions::VarCoeffMap &dirForcing, Array< OneD, NekDouble > &outarray)
void AddHDGHelmholtzTraceTerms (const NekDouble tau, const Array< OneD, const NekDouble > &inarray, Array< OneD, ExpansionSharedPtr > &EdgeExp, const StdRegions::VarCoeffMap &dirForcing, Array< OneD, NekDouble > &outarray)
Expansion3DSharedPtr GetLeftAdjacentElementExp () const
Expansion3DSharedPtr GetRightAdjacentElementExp () const
int GetLeftAdjacentElementFace () const
int GetRightAdjacentElementFace () const
void SetAdjacentElementExp (int face, Expansion3DSharedPtr &f)
SpatialDomains::Geometry2DSharedPtr GetGeom2D () 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
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)
 Integrates the specified function over the domain.
virtual NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
virtual void v_PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 Calculate the derivative of the physical points.
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.
virtual void v_PhysDirectionalDeriv (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &direction, Array< OneD, NekDouble > &out)
 Physical derivative along a direction vector.
virtual void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Transform a given function from physical quadrature space to coefficient space.
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)
 Calculate the inner product of inarray with respect to the basis B=base0*base1 and put into outarray.
virtual void v_IProductWRTDerivBase (const int dir, 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)
virtual void v_IProductWRTBase_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_IProductWRTDerivBase_MatOp (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
virtual void v_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 NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &coord, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain.
virtual void v_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.
virtual void v_GetEdgePhysVals (const int edge, const StdRegions::StdExpansionSharedPtr &EdgeExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, 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 void v_GetEdgeInterpVals (const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_GetEdgeQFactors (const int edge, Array< OneD, NekDouble > &outarray)
virtual void v_ComputeEdgeNormal (const int edge)
virtual const
SpatialDomains::GeomFactorsSharedPtr
v_GetMetricInfo () const
virtual int v_GetCoordim ()
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 the same expansion type.
virtual StdRegions::Orientation v_GetEorient (int edge)
virtual StdRegions::Orientation v_GetCartesianEorient (int edge)
virtual const
LibUtilities::BasisSharedPtr
v_GetBasis (int dir) const
virtual int v_GetNumPoints (const int dir) const
virtual DNekMatSharedPtr v_GenMatrix (const StdRegions::StdMatrixKey &mkey)
virtual DNekMatSharedPtr v_CreateStdMatrix (const StdRegions::StdMatrixKey &mkey)
DNekScalMatSharedPtr CreateMatrix (const MatrixKey &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)
virtual void v_MassMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
virtual void v_LaplacianMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
virtual void v_LaplacianMatrixOp (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
virtual void v_WeakDerivMatrixOp (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
virtual void v_WeakDirectionalDerivMatrixOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
virtual void v_MassLevelCurvatureMatrixOp (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 void v_GeneralMatrixOp_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
virtual void v_LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
virtual void v_ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_ComputeLaplacianMetric ()
- Protected Member Functions inherited from Nektar::StdRegions::StdQuadExp
virtual void v_StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
virtual void v_StdPhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_BwdTrans_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1)
virtual void v_IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1)
virtual void v_LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
virtual void v_FillMode (const int mode, Array< OneD, NekDouble > &array)
 Fill outarray with mode mode of expansion.
virtual int v_GetNverts () const
virtual int v_GetNedges () const
virtual int v_GetEdgeNcoeffs (const int i) const
virtual int v_GetEdgeNumPoints (const int i) const
virtual int v_NumBndryCoeffs () const
virtual int v_NumDGBndryCoeffs () const
virtual LibUtilities::BasisType v_GetEdgeBasisType (const int i) const
virtual int v_DetCartesianDirOfEdge (const int edge)
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
virtual const
LibUtilities::BasisKey 
v_DetEdgeBasisKey (const int i) const
virtual LibUtilities::ShapeType v_DetShapeType () const
virtual bool v_IsBoundaryInteriorExpansion ()
virtual void v_GetBoundaryMap (Array< OneD, unsigned int > &outarray)
virtual void v_GetInteriorMap (Array< OneD, unsigned int > &outarray)
virtual int v_GetVertexMap (int localVertexId, bool useCoeffPacking=false)
virtual void v_GetEdgeInteriorMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
void v_GetEdgeToElementMap (const int eid, const Orientation edgeOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
virtual void v_SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_GetSimplexEquiSpacedConnectivity (Array< OneD, int > &conn, bool standard=true)
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion2D
virtual NekDouble v_PhysEvaluate (const Array< OneD, DNekMatSharedPtr > &I, const Array< OneD, const NekDouble > &physvals)
virtual void v_LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
virtual void v_HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion
DNekMatSharedPtr CreateStdMatrix (const StdMatrixKey &mkey)
DNekBlkMatSharedPtr CreateStdStaticCondMatrix (const StdMatrixKey &mkey)
 Create the static condensation of a matrix when using a boundary interior decomposition.
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.
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 IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void GeneralMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void MassMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LaplacianMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
void LaplacianMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LaplacianMatrixOp_MatFree (const int k1, const int k2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void WeakDerivMatrixOp_MatFree (const int i, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void WeakDirectionalDerivMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void MassLevelCurvatureMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LinearAdvectionDiffusionReactionMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey, bool addDiffusionTerm=true)
void HelmholtzMatrixOp_MatFree (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void HelmholtzMatrixOp_MatFree_GenericImpl (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
virtual void v_SetCoeffsToOrientation (StdRegions::Orientation dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_SetCoeffsToOrientation (Array< OneD, NekDouble > &coeffs, StdRegions::Orientation dir)
- Protected Member Functions inherited from Nektar::LocalRegions::Expansion2D
virtual void v_DGDeriv (const int dir, const Array< OneD, const NekDouble > &incoeffs, Array< OneD, ExpansionSharedPtr > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &edgeCoeffs, Array< OneD, NekDouble > &out_d)
virtual void v_AddEdgeNormBoundaryInt (const int edge, const ExpansionSharedPtr &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 ExpansionSharedPtr &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
virtual void v_AddRobinMassMatrix (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
virtual void v_AddRobinEdgeContribution (const int edgeid, const Array< OneD, const NekDouble > &primCoeffs, Array< OneD, NekDouble > &coeffs)
virtual DNekMatSharedPtr v_BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd)
void GetPhysEdgeVarCoeffsFromElement (const int edge, ExpansionSharedPtr &EdgeExp, const Array< OneD, const NekDouble > &varcoeff, Array< OneD, NekDouble > &outarray)
Array< OneD, unsigned int > v_GetEdgeInverseBoundaryMap (int eid)
virtual void v_NegateEdgeNormal (const int edge)
virtual bool v_EdgeNormalNegated (const int edge)
virtual void v_SetUpPhysNormals (const int edge)
const StdRegions::NormalVectorv_GetEdgeNormal (const int edge) const
const StdRegions::NormalVectorv_GetSurfaceNormal (const int id) const
- 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 DNekMatSharedPtr v_BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
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)

Private Member Functions

 QuadExp ()

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::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::Expansion2D
std::vector< ExpansionWeakPtrm_edgeExp
std::vector< bool > m_requireNeg
std::map< int,
StdRegions::NormalVector
m_edgeNormals
std::map< int, bool > m_negatedNormals
Expansion3DWeakPtr m_elementLeft
Expansion3DWeakPtr m_elementRight
int m_elementFaceLeft
int m_elementFaceRight

Detailed Description

Definition at line 52 of file QuadExp.h.

Constructor & Destructor Documentation

Nektar::LocalRegions::QuadExp::QuadExp ( const LibUtilities::BasisKey Ba,
const LibUtilities::BasisKey Bb,
const SpatialDomains::QuadGeomSharedPtr geom 
)

Constructor using BasisKey class for quadrature points and order definition.

Definition at line 48 of file QuadExp.cpp.

:
StdExpansion (Ba.GetNumModes()*Bb.GetNumModes(),2,Ba,Bb),
StdExpansion2D(Ba.GetNumModes()*Bb.GetNumModes(),Ba,Bb),
StdQuadExp (Ba,Bb),
Expansion (geom),
Expansion2D (geom),
boost::bind(&QuadExp::CreateMatrix, this, _1),
std::string("QuadExpMatrix")),
boost::bind(&QuadExp::CreateStaticCondMatrix, this, _1),
std::string("QuadExpStaticCondMatrix"))
{
}
Nektar::LocalRegions::QuadExp::QuadExp ( const QuadExp T)

Definition at line 66 of file QuadExp.cpp.

:
Expansion (T),
m_matrixManager(T.m_matrixManager),
m_staticCondMatrixManager(T.m_staticCondMatrixManager)
{
}
Nektar::LocalRegions::QuadExp::~QuadExp ( )
virtual

Definition at line 78 of file QuadExp.cpp.

{
}
Nektar::LocalRegions::QuadExp::QuadExp ( )
private

Member Function Documentation

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

Definition at line 1533 of file QuadExp.cpp.

References ASSERTL2, Nektar::LocalRegions::Expansion::BuildVertexMatrix(), Nektar::StdRegions::StdMatrixKey::ConstFactorExists(), Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::eFactorGaussEdge, Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::eHelmholtz, Nektar::StdRegions::eHybridDGHelmholtz, Nektar::StdRegions::eInterpGauss, Nektar::StdRegions::eInvHybridDGHelmholtz, Nektar::StdRegions::eInvLaplacianWithUnityMean, Nektar::StdRegions::eInvMass, Nektar::StdRegions::eIProductWRTBase, Nektar::StdRegions::eIProductWRTDerivBase0, Nektar::StdRegions::eIProductWRTDerivBase1, Nektar::StdRegions::eIProductWRTDerivBase2, Nektar::StdRegions::eLaplacian, Nektar::StdRegions::eLaplacian00, Nektar::StdRegions::eLaplacian01, Nektar::StdRegions::eLaplacian11, Nektar::StdRegions::eMass, Nektar::SpatialDomains::eNoGeomType, Nektar::StdRegions::ePreconLinearSpace, Nektar::StdRegions::eWeakDeriv0, Nektar::StdRegions::eWeakDeriv1, Nektar::StdRegions::eWeakDeriv2, Nektar::StdRegions::StdExpansion::GenMatrix(), Nektar::StdRegions::StdMatrixKey::GetConstFactor(), Nektar::StdRegions::StdMatrixKey::GetConstFactors(), Nektar::StdRegions::StdExpansion::GetLocStaticCondMatrix(), 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, m_matrixManager, Nektar::LocalRegions::Expansion::m_metricinfo, and Nektar::Transpose().

{
"Geometric information is not set up");
switch (mkey.GetMatrixType())
{
{
if ((m_metricinfo->GetGtype() ==
SpatialDomains::eDeformed) || (mkey.GetNVarCoeff()))
{
NekDouble one = 1.0;
returnval = MemoryManager<DNekScalMat>::
AllocateSharedPtr(one,mat);
}
else
{
NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
returnval = MemoryManager<DNekScalMat>::
AllocateSharedPtr(jac,mat);
}
}
break;
{
if ((m_metricinfo->GetGtype() ==
SpatialDomains::eDeformed) || (mkey.GetNVarCoeff()))
{
NekDouble one = 1.0;
StdRegions::StdMatrixKey masskey(
DNekMatSharedPtr mat = GenMatrix(masskey);
mat->Invert();
returnval = MemoryManager<DNekScalMat>::
AllocateSharedPtr(one,mat);
}
else
{
NekDouble fac = 1.0/(m_metricinfo->GetJac(ptsKeys))[0];
returnval = MemoryManager<DNekScalMat>::
AllocateSharedPtr(fac,mat);
}
}
break;
{
|| (mkey.GetNVarCoeff()))
{
NekDouble one = 1.0;
returnval = MemoryManager<DNekScalMat>::
AllocateSharedPtr(one,mat);
}
else
{
NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
Array<TwoD, const NekDouble> df =
m_metricinfo->GetDerivFactors(ptsKeys);
int dir = 0;
switch(mkey.GetMatrixType())
{
dir = 0;
break;
dir = 1;
break;
dir = 2;
break;
default:
break;
}
MatrixKey deriv0key(StdRegions::eWeakDeriv0,
mkey.GetShapeType(), *this);
MatrixKey deriv1key(StdRegions::eWeakDeriv1,
mkey.GetShapeType(), *this);
DNekMat &deriv0 = *GetStdMatrix(deriv0key);
DNekMat &deriv1 = *GetStdMatrix(deriv1key);
int rows = deriv0.GetRows();
int cols = deriv1.GetColumns();
DNekMatSharedPtr WeakDeriv = MemoryManager<DNekMat>::
AllocateSharedPtr(rows,cols);
(*WeakDeriv) = df[2*dir][0]*deriv0 +
df[2*dir+1][0]*deriv1;
returnval = MemoryManager<DNekScalMat>::
AllocateSharedPtr(jac,WeakDeriv);
}
}
break;
{
if( (m_metricinfo->GetGtype() ==
SpatialDomains::eDeformed) || (mkey.GetNVarCoeff() > 0)
|| (mkey.ConstFactorExists
{
NekDouble one = 1.0;
returnval = MemoryManager<DNekScalMat>::
AllocateSharedPtr(one,mat);
}
else
{
MatrixKey lap00key(StdRegions::eLaplacian00,
mkey.GetShapeType(), *this);
MatrixKey lap01key(StdRegions::eLaplacian01,
mkey.GetShapeType(), *this);
MatrixKey lap11key(StdRegions::eLaplacian11,
mkey.GetShapeType(), *this);
DNekMat &lap00 = *GetStdMatrix(lap00key);
DNekMat &lap01 = *GetStdMatrix(lap01key);
DNekMat &lap11 = *GetStdMatrix(lap11key);
NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
Array<TwoD, const NekDouble>
gmat = m_metricinfo->GetGmat(ptsKeys);
int rows = lap00.GetRows();
int cols = lap00.GetColumns();
MemoryManager<DNekMat>::AllocateSharedPtr(rows,cols);
(*lap) = gmat[0][0] * lap00 +
gmat[1][0] * (lap01 + Transpose(lap01)) +
gmat[3][0] * lap11;
returnval = MemoryManager<DNekScalMat>::
AllocateSharedPtr(jac,lap);
}
}
break;
{
returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(1.0,mat);
}
break;
{
NekDouble lambda =
mkey.GetConstFactor(StdRegions::eFactorLambda);
MatrixKey masskey(mkey, StdRegions::eMass);
DNekScalMat &MassMat = *(this->m_matrixManager[masskey]);
MatrixKey lapkey(mkey, StdRegions::eLaplacian);
DNekScalMat &LapMat = *(this->m_matrixManager[lapkey]);
int rows = LapMat.GetRows();
int cols = LapMat.GetColumns();
DNekMatSharedPtr helm = MemoryManager<DNekMat>::
AllocateSharedPtr(rows,cols);
NekDouble one = 1.0;
(*helm) = LapMat + lambda*MassMat;
returnval =
MemoryManager<DNekScalMat>::AllocateSharedPtr(one,helm);
}
break;
{
{
NekDouble one = 1.0;
returnval = MemoryManager<DNekScalMat>::
AllocateSharedPtr(one,mat);
}
else
{
NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
returnval = MemoryManager<DNekScalMat>::
AllocateSharedPtr(jac,mat);
}
}
break;
{
{
NekDouble one = 1.0;
returnval = MemoryManager<DNekScalMat>::
AllocateSharedPtr(one,mat);
}
else
{
NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
const Array<TwoD, const NekDouble>& df =
m_metricinfo->GetDerivFactors(ptsKeys);
int dir = 0;
switch(mkey.GetMatrixType())
{
dir = 0;
break;
dir = 1;
break;
dir = 2;
break;
default:
break;
}
MatrixKey iProdDeriv0Key(
mkey.GetShapeType(), *this);
MatrixKey iProdDeriv1Key(
mkey.GetShapeType(), *this);
DNekMat &stdiprod0 = *GetStdMatrix(iProdDeriv0Key);
DNekMat &stdiprod1 = *GetStdMatrix(iProdDeriv0Key);
int rows = stdiprod0.GetRows();
int cols = stdiprod1.GetColumns();
DNekMatSharedPtr mat = MemoryManager<DNekMat>::
AllocateSharedPtr(rows,cols);
(*mat) = df[2*dir][0]*stdiprod0 +
df[2*dir+1][0]*stdiprod1;
returnval = MemoryManager<DNekScalMat>::
AllocateSharedPtr(jac,mat);
}
}
break;
{
NekDouble one = 1.0;
DetShapeType(), *this,
mkey.GetConstFactors(), mkey.GetVarCoeffs());
mat->Invert();
returnval =
MemoryManager<DNekScalMat>::AllocateSharedPtr(one,mat);
}
break;
{
Array<OneD, NekDouble> coords(1, 0.0);
StdRegions::ConstFactorMap factors = mkey.GetConstFactors();
int edge = (int)factors[StdRegions::eFactorGaussEdge];
coords[0] = (edge == 0 || edge == 3) ? -1.0 : 1.0;
m_Ix = m_base[(edge + 1) % 2]->GetI(coords);
returnval =
MemoryManager<DNekScalMat>::AllocateSharedPtr(1.0,m_Ix);
}
break;
{
NekDouble one = 1.0;
MatrixKey helmkey(
StdRegions::eHelmholtz, mkey.GetShapeType(), *this,
mkey.GetConstFactors(), mkey.GetVarCoeffs());
DNekScalBlkMatSharedPtr helmStatCond =
DNekScalMatSharedPtr A =helmStatCond->GetBlock(0,0);
returnval =
MemoryManager<DNekScalMat>::AllocateSharedPtr(one, R);
}
break;
default:
{
NekDouble one = 1.0;
returnval =
MemoryManager<DNekScalMat>::AllocateSharedPtr(one, mat);
}
break;
}
return returnval;
}
DNekScalBlkMatSharedPtr Nektar::LocalRegions::QuadExp::CreateStaticCondMatrix ( const MatrixKey mkey)
protected

Definition at line 1847 of file QuadExp.cpp.

References ASSERTL2, Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::eMass, Nektar::SpatialDomains::eNoGeomType, Nektar::StdRegions::StdExpansion::GetBoundaryMap(), Nektar::StdRegions::StdExpansion::GetInteriorMap(), Nektar::LocalRegions::Expansion::GetLocMatrix(), Nektar::StdRegions::StdMatrixKey::GetMatrixType(), Nektar::StdRegions::StdMatrixKey::GetNVarCoeff(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::GetStdStaticCondMatrix(), Nektar::LocalRegions::Expansion::m_metricinfo, Nektar::StdRegions::StdExpansion::m_ncoeffs, and Nektar::StdRegions::StdExpansion::NumBndryCoeffs().

{
ASSERTL2(m_metricinfo->GetGtype()
"Geometric information is not set up");
// set up block matrix system
unsigned int nbdry = NumBndryCoeffs();
unsigned int nint = (unsigned int)(m_ncoeffs - nbdry);
unsigned int exp_size[] = {nbdry,nint};
unsigned int nblks = 2;
returnval = MemoryManager<DNekScalBlkMat>::
AllocateSharedPtr(nblks,nblks,exp_size,exp_size);
//Really need a constructor which takes Arrays
NekDouble factor = 1.0;
switch (mkey.GetMatrixType())
{
// this can only use stdregions statically condensed system
// for mass matrix
||(mkey.GetNVarCoeff()))
{
factor = 1.0;
goto UseLocRegionsMatrix;
}
else
{
factor = (m_metricinfo->GetJac(GetPointsKeys()))[0];
goto UseStdRegionsMatrix;
}
break;
default: // use Deformed case for both
// regular and deformed geometries
factor = 1.0;
goto UseLocRegionsMatrix;
break;
UseStdRegionsMatrix:
{
NekDouble invfactor = 1.0/factor;
NekDouble one = 1.0;
returnval->SetBlock(0,0,Atmp = MemoryManager<DNekScalMat>::
AllocateSharedPtr(factor,Asubmat = mat->GetBlock(0,0)));
returnval->SetBlock(0,1,Atmp = MemoryManager<DNekScalMat>::
AllocateSharedPtr(one,Asubmat = mat->GetBlock(0,1)));
returnval->SetBlock(1,0,Atmp = MemoryManager<DNekScalMat>::
AllocateSharedPtr(factor,Asubmat = mat->GetBlock(1,0)));
returnval->SetBlock(1,1,Atmp = MemoryManager<DNekScalMat>::
AllocateSharedPtr(invfactor,Asubmat = mat->GetBlock(1,1)));
}
break;
UseLocRegionsMatrix:
{
int i,j;
NekDouble invfactor = 1.0/factor;
NekDouble one = 1.0;
DNekScalMat &mat = *GetLocMatrix(mkey);
DNekMatSharedPtr A = MemoryManager<DNekMat>::
AllocateSharedPtr(nbdry,nbdry);
DNekMatSharedPtr B = MemoryManager<DNekMat>::
AllocateSharedPtr(nbdry,nint);
DNekMatSharedPtr C = MemoryManager<DNekMat>::
AllocateSharedPtr(nint,nbdry);
DNekMatSharedPtr D = MemoryManager<DNekMat>::
AllocateSharedPtr(nint,nint);
Array<OneD,unsigned int> bmap(nbdry);
Array<OneD,unsigned int> imap(nint);
for (i = 0; i < nbdry; ++i)
{
for(j = 0; j < nbdry; ++j)
{
(*A)(i,j) = mat(bmap[i],bmap[j]);
}
for(j = 0; j < nint; ++j)
{
(*B)(i,j) = mat(bmap[i],imap[j]);
}
}
for (i = 0; i < nint; ++i)
{
for(j = 0; j < nbdry; ++j)
{
(*C)(i,j) = mat(imap[i],bmap[j]);
}
for(j = 0; j < nint; ++j)
{
(*D)(i,j) = mat(imap[i],imap[j]);
}
}
// Calculate static condensed system
if(nint)
{
D->Invert();
(*B) = (*B)*(*D);
(*A) = (*A) - (*B)*(*C);
}
returnval->SetBlock(0,0,Atmp = MemoryManager<DNekScalMat>::
AllocateSharedPtr(factor, A));
returnval->SetBlock(0,1,Atmp = MemoryManager<DNekScalMat>::
AllocateSharedPtr(one, B));
returnval->SetBlock(1,0,Atmp = MemoryManager<DNekScalMat>::
AllocateSharedPtr(factor, C));
returnval->SetBlock(1,1,Atmp = MemoryManager<DNekScalMat>::
AllocateSharedPtr(invfactor, D));
}
}
return returnval;
}
void Nektar::LocalRegions::QuadExp::v_ComputeEdgeNormal ( const int  edge)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1116 of file QuadExp.cpp.

References ASSERTL0, Nektar::StdRegions::eBackwards, Nektar::SpatialDomains::eDeformed, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::SpatialDomains::eMovingRegular, Nektar::SpatialDomains::eRegular, Vmath::Fill(), Nektar::StdRegions::StdExpansion::GetCoordim(), Nektar::StdRegions::StdExpansion::GetEorient(), Nektar::LocalRegions::Expansion::GetGeom(), Nektar::LibUtilities::PointsKey::GetNumPoints(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::LibUtilities::Interp1D(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion2D::m_edgeNormals, Vmath::Reverse(), Vmath::Sdiv(), Vmath::Smul(), v_GetEdgeInterpVals(), Vmath::Vmul(), Vmath::Vsqrt(), Vmath::Vvtvp(), and Vmath::Zero().

{
int i;
GetGeom()->GetMetricInfo();
SpatialDomains::GeomType type = geomFactors->GetGtype();
const Array<TwoD, const NekDouble> & df = geomFactors->GetDerivFactors(ptsKeys);
const Array<OneD, const NekDouble> & jac = geomFactors->GetJac(ptsKeys);
int nqe = m_base[0]->GetNumPoints();
int vCoordDim = GetCoordim();
m_edgeNormals[edge] = Array<OneD, Array<OneD, NekDouble> >
(vCoordDim);
Array<OneD, Array<OneD, NekDouble> > &normal = m_edgeNormals[edge];
for (i = 0; i < vCoordDim; ++i)
{
normal[i] = Array<OneD, NekDouble>(nqe);
}
// Regular geometry case
if ((type == SpatialDomains::eRegular)||
{
NekDouble fac;
// Set up normals
switch (edge)
{
case 0:
for (i = 0; i < vCoordDim; ++i)
{
Vmath::Fill(nqe, -df[2*i+1][0], normal[i], 1);
}
break;
case 1:
for (i = 0; i < vCoordDim; ++i)
{
Vmath::Fill(nqe, df[2*i][0], normal[i], 1);
}
break;
case 2:
for (i = 0; i < vCoordDim; ++i)
{
Vmath::Fill(nqe, df[2*i+1][0], normal[i], 1);
}
break;
case 3:
for (i = 0; i < vCoordDim; ++i)
{
Vmath::Fill(nqe, -df[2*i][0], normal[i], 1);
}
break;
default:
ASSERTL0(false, "edge is out of range (edge < 4)");
}
// normalise
fac = 0.0;
for (i =0 ; i < vCoordDim; ++i)
{
fac += normal[i][0]*normal[i][0];
}
fac = 1.0/sqrt(fac);
for (i = 0; i < vCoordDim; ++i)
{
Vmath::Smul(nqe, fac, normal[i], 1,normal[i], 1);
}
}
else // Set up deformed normals
{
int j;
int nquad0 = ptsKeys[0].GetNumPoints();
int nquad1 = ptsKeys[1].GetNumPoints();
LibUtilities::PointsKey from_key;
Array<OneD,NekDouble> normals(vCoordDim*max(nquad0,nquad1),0.0);
Array<OneD,NekDouble> edgejac(vCoordDim*max(nquad0,nquad1),0.0);
// Extract Jacobian along edges and recover local
// derivates (dx/dr) for polynomial interpolation by
// multiplying m_gmat by jacobian
// Implementation for all the basis except Gauss points
if (m_base[0]->GetPointsType() !=
&& m_base[1]->GetPointsType() !=
{
switch (edge)
{
case 0:
for (j = 0; j < nquad0; ++j)
{
edgejac[j] = jac[j];
for (i = 0; i < vCoordDim; ++i)
{
normals[i*nquad0+j] =
-df[2*i+1][j]*edgejac[j];
}
}
from_key = ptsKeys[0];
break;
case 1:
for (j = 0; j < nquad1; ++j)
{
edgejac[j] = jac[nquad0*j+nquad0-1];
for (i = 0; i < vCoordDim; ++i)
{
normals[i*nquad1+j] =
df[2*i][nquad0*j + nquad0-1]
*edgejac[j];
}
}
from_key = ptsKeys[1];
break;
case 2:
for (j = 0; j < nquad0; ++j)
{
edgejac[j] = jac[nquad0*(nquad1-1)+j];
for (i = 0; i < vCoordDim; ++i)
{
normals[i*nquad0+j] =
(df[2*i+1][nquad0*(nquad1-1)+j])
*edgejac[j];
}
}
from_key = ptsKeys[0];
break;
case 3:
for (j = 0; j < nquad1; ++j)
{
edgejac[j] = jac[nquad0*j];
for (i = 0; i < vCoordDim; ++i)
{
normals[i*nquad1+j] =
-df[2*i][nquad0*j]*edgejac[j];
}
}
from_key = ptsKeys[1];
break;
default:
ASSERTL0(false,"edge is out of range (edge < 3)");
}
}
else
{
int nqtot = nquad0 * nquad1;
Array<OneD, NekDouble> tmp_gmat(nqtot, 0.0);
Array<OneD, NekDouble> tmp_gmat_edge(nqe, 0.0);
switch (edge)
{
case 0:
for (j = 0; j < nquad0; ++j)
{
for (i = 0; i < vCoordDim; ++i)
{
Vmath::Vmul(nqtot,
&(df[2*i+1][0]), 1,
&jac[0], 1,
&(tmp_gmat[0]), 1);
edge, tmp_gmat, tmp_gmat_edge);
normals[i*nquad0+j] = -tmp_gmat_edge[j];
}
}
from_key = ptsKeys[0];
break;
case 1:
for (j = 0; j < nquad1; ++j)
{
for (i = 0; i < vCoordDim; ++i)
{
Vmath::Vmul(nqtot,
&(df[2*i][0]), 1,
&jac[0], 1,
&(tmp_gmat[0]), 1);
edge, tmp_gmat, tmp_gmat_edge);
normals[i*nquad1+j] = tmp_gmat_edge[j];
}
}
from_key = ptsKeys[1];
break;
case 2:
for (j = 0; j < nquad0; ++j)
{
for (i = 0; i < vCoordDim; ++i)
{
Vmath::Vmul(nqtot,
&(df[2*i+1][0]), 1,
&jac[0], 1,
&(tmp_gmat[0]), 1);
edge, tmp_gmat, tmp_gmat_edge);
normals[i*nquad0+j] = tmp_gmat_edge[j];
}
}
from_key = ptsKeys[0];
break;
case 3:
for (j = 0; j < nquad1; ++j)
{
for (i = 0; i < vCoordDim; ++i)
{
Vmath::Vmul(nqtot,
&(df[2*i][0]), 1,
&jac[0], 1,
&(tmp_gmat[0]) ,1);
edge, tmp_gmat, tmp_gmat_edge);
normals[i*nquad1+j] = -tmp_gmat_edge[j];
}
}
from_key = ptsKeys[1];
break;
default:
ASSERTL0(false,"edge is out of range (edge < 3)");
}
}
int nq = from_key.GetNumPoints();
Array<OneD,NekDouble> work(nqe,0.0);
// interpolate Jacobian and invert
from_key,jac, m_base[0]->GetPointsKey(), work);
Vmath::Sdiv(nq,1.0,&work[0],1,&work[0],1);
// interpolate
for (i = 0; i < GetCoordim(); ++i)
{
from_key,&normals[i*nq],
m_base[0]->GetPointsKey(),
&normal[i][0]);
Vmath::Vmul(nqe, work, 1, normal[i], 1, normal[i], 1);
}
//normalise normal vectors
Vmath::Zero(nqe,work,1);
for (i = 0; i < GetCoordim(); ++i)
{
normal[i], 1,
normal[i],1 ,
work, 1,
work, 1);
}
Vmath::Vsqrt(nqe,work,1,work,1);
Vmath::Sdiv(nqe,1.0,work,1,work,1);
for (i = 0; i < GetCoordim(); ++i)
{
Vmath::Vmul(nqe, normal[i], 1, work, 1, normal[i], 1);
}
// Reverse direction so that points are in
// anticlockwise direction if edge >=2
if (edge >= 2)
{
for (i = 0; i < GetCoordim(); ++i)
{
Vmath::Reverse(nqe, normal[i], 1, normal[i], 1);
}
}
}
{
for (i = 0; i < vCoordDim; ++i)
{
if (geomFactors->GetGtype() == SpatialDomains::eDeformed)
{
Vmath::Reverse(nqe, normal[i], 1, normal[i],1);
}
}
}
}
void Nektar::LocalRegions::QuadExp::v_ComputeLaplacianMetric ( )
protectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 2191 of file QuadExp.cpp.

References Nektar::LocalRegions::Expansion::ComputeQuadratureMetric(), Nektar::SpatialDomains::eDeformed, Vmath::Fill(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::GetTotPoints(), Nektar::LocalRegions::Expansion::m_metricinfo, Nektar::LocalRegions::Expansion::m_metrics, Nektar::LocalRegions::MetricLaplacian00, Nektar::LocalRegions::MetricLaplacian01, Nektar::LocalRegions::MetricLaplacian02, Nektar::LocalRegions::MetricLaplacian11, Nektar::LocalRegions::MetricLaplacian12, Nektar::LocalRegions::MetricLaplacian22, Nektar::LocalRegions::MetricQuadrature, Nektar::StdRegions::StdExpansion::MultiplyByQuadratureMetric(), and Vmath::Vcopy().

{
if (m_metrics.count(MetricQuadrature) == 0)
{
}
const SpatialDomains::GeomType type = m_metricinfo->GetGtype();
const unsigned int nqtot = GetTotPoints();
const unsigned int dim = 2;
};
const Array<TwoD, const NekDouble> gmat =
for (unsigned int i = 0; i < dim; ++i)
{
for (unsigned int j = i; j < dim; ++j)
{
m_metrics[m[i][j]] = Array<OneD, NekDouble>(nqtot);
{
Vmath::Vcopy(nqtot, &gmat[i*dim+j][0], 1,
&m_metrics[m[i][j]][0], 1);
}
else
{
Vmath::Fill(nqtot, gmat[i*dim+j][0],
&m_metrics[m[i][j]][0], 1);
}
m_metrics[m[i][j]]);
}
}
}
DNekMatSharedPtr Nektar::LocalRegions::QuadExp::v_CreateStdMatrix ( const StdRegions::StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdQuadExp.

Definition at line 1522 of file QuadExp.cpp.

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

{
LibUtilities::BasisKey bkey0 = m_base[0]->GetBasisKey();
LibUtilities::BasisKey bkey1 = m_base[1]->GetBasisKey();
MemoryManager<StdQuadExp>::AllocateSharedPtr(bkey0,bkey1);
return tmp->GetStdMatrix(mkey);
}
void Nektar::LocalRegions::QuadExp::v_DropLocStaticCondMatrix ( const MatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1989 of file QuadExp.cpp.

References Nektar::LibUtilities::NekManager< KeyType, ValueT, opLessCreator >::DeleteObject(), and m_staticCondMatrixManager.

void Nektar::LocalRegions::QuadExp::v_ExtractDataToCoeffs ( const NekDouble data,
const std::vector< unsigned int > &  nummodes,
const int  nmode_offset,
NekDouble coeffs 
)
protectedvirtual

Unpack data from input file assuming it comes from the same expansion type.

See Also
StdExpansion::ExtractDataToCoeffs

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1410 of file QuadExp.cpp.

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

{
int data_order0 = nummodes[mode_offset];
int fillorder0 = std::min(m_base[0]->GetNumModes(),data_order0);
int data_order1 = nummodes[mode_offset + 1];
int order1 = m_base[1]->GetNumModes();
int fillorder1 = min(order1,data_order1);
switch (m_base[0]->GetBasisType())
{
{
int i;
int cnt = 0;
int cnt1 = 0;
"Extraction routine not set up for this basis");
for (i = 0; i < fillorder0; ++i)
{
Vmath::Vcopy(fillorder1, data + cnt, 1, coeffs +cnt1, 1);
cnt += data_order1;
cnt1 += order1;
}
}
break;
{
// Assume that input is also Gll_Lagrange but no way to check;
LibUtilities::PointsKey
LibUtilities::PointsKey
LibUtilities::Interp2D(p0, p1, data,
m_base[0]->GetPointsKey(),
m_base[1]->GetPointsKey(),
coeffs);
}
break;
{
// Assume that input is also Gll_Lagrange but no way to check;
LibUtilities::PointsKey
LibUtilities::PointsKey
LibUtilities::Interp2D(p0, p1, data,
m_base[0]->GetPointsKey(),
m_base[1]->GetPointsKey(),
coeffs);
}
break;
default:
ASSERTL0(false,
"basis is either not set up or not hierarchicial");
}
}
void Nektar::LocalRegions::QuadExp::v_FwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Transform a given function from physical quadrature space to coefficient space.

See Also
StdExpansion::FwdTrans

Reimplemented from Nektar::StdRegions::StdQuadExp.

Definition at line 260 of file QuadExp.cpp.

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

{
if ((m_base[0]->Collocation())&&(m_base[1]->Collocation()))
{
Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
}
else
{
IProductWRTBase(inarray,outarray);
// get Mass matrix inverse
MatrixKey masskey(StdRegions::eInvMass,
DetShapeType(),*this);
// copy inarray in case inarray == outarray
NekVector<NekDouble> in(m_ncoeffs,outarray,eCopy);
NekVector<NekDouble> out(m_ncoeffs,outarray,eWrapper);
out = (*matsys)*in;
}
}
void Nektar::LocalRegions::QuadExp::v_FwdTrans_BndConstrained ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdQuadExp.

Definition at line 286 of file QuadExp.cpp.

References Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::StdRegions::eBackwards, Nektar::StdRegions::eMass, GetEdge(), Nektar::StdRegions::StdExpansion::GetEdgeToElementMap(), Nektar::StdRegions::StdExpansion::GetEorient(), Nektar::LocalRegions::Expansion2D::GetGeom2D(), Nektar::StdRegions::StdExpansion::GetInteriorMap(), Nektar::StdRegions::StdExpansion::IProductWRTBase(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, m_staticCondMatrixManager, Nektar::StdRegions::StdExpansion::MassMatrixOp(), Nektar::StdRegions::StdExpansion::NumBndryCoeffs(), sign, Vmath::Vcopy(), and Vmath::Vsub().

{
if ((m_base[0]->Collocation())&&(m_base[1]->Collocation()))
{
Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
}
else
{
int i,j;
int npoints[2] = {m_base[0]->GetNumPoints(),
m_base[1]->GetNumPoints()};
int nmodes[2] = {m_base[0]->GetNumModes(),
m_base[1]->GetNumModes()};
fill(outarray.get(), outarray.get()+m_ncoeffs, 0.0 );
Array<OneD, NekDouble> physEdge[4];
Array<OneD, NekDouble> coeffEdge[4];
for (i = 0; i < 4; i++)
{
physEdge[i] = Array<OneD, NekDouble>(npoints[i%2]);
coeffEdge[i] = Array<OneD, NekDouble>(nmodes[i%2]);
orient[i] = GetEorient(i);
}
for (i = 0; i < npoints[0]; i++)
{
physEdge[0][i] = inarray[i];
physEdge[2][i] = inarray[npoints[0]*npoints[1]-1-i];
}
for (i = 0; i < npoints[1]; i++)
{
physEdge[1][i] =
inarray[npoints[0]-1+i*npoints[0]];
physEdge[3][i] =
inarray[(npoints[1]-1)*npoints[0]-i*npoints[0]];
}
for (i = 0; i < 4; i++)
{
if ( orient[i] == StdRegions::eBackwards )
{
reverse((physEdge[i]).get(),
(physEdge[i]).get() + npoints[i%2] );
}
}
SegExpSharedPtr segexp[4];
for (i = 0; i < 4; i++)
{
segexp[i] = MemoryManager<LocalRegions::SegExp>::
AllocateSharedPtr(
m_base[i%2]->GetBasisKey(),GetGeom2D()->GetEdge(i));
}
Array<OneD, unsigned int> mapArray;
Array<OneD, int> signArray;
for (i = 0; i < 4; i++)
{
segexp[i%2]->FwdTrans_BndConstrained(
physEdge[i],coeffEdge[i]);
GetEdgeToElementMap(i,orient[i],mapArray,signArray);
for (j=0; j < nmodes[i%2]; j++)
{
sign = (NekDouble) signArray[j];
outarray[ mapArray[j] ] = sign * coeffEdge[i][j];
}
}
int nBoundaryDofs = NumBndryCoeffs();
int nInteriorDofs = m_ncoeffs - nBoundaryDofs;
if (nInteriorDofs > 0) {
Array<OneD, NekDouble> tmp0(m_ncoeffs);
Array<OneD, NekDouble> tmp1(m_ncoeffs);
StdRegions::StdMatrixKey
stdmasskey(StdRegions::eMass,DetShapeType(),*this);
MassMatrixOp(outarray,tmp0,stdmasskey);
IProductWRTBase(inarray,tmp1);
Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
// get Mass matrix inverse (only of interior DOF)
// use block (1,1) of the static condensed system
// note: this block alreay contains the inverse matrix
MatrixKey
masskey(StdRegions::eMass,DetShapeType(),*this);
(m_staticCondMatrixManager[masskey])->GetBlock(1,1);
Array<OneD, NekDouble> rhs(nInteriorDofs);
Array<OneD, NekDouble> result(nInteriorDofs);
GetInteriorMap(mapArray);
for (i = 0; i < nInteriorDofs; i++)
{
rhs[i] = tmp1[ mapArray[i] ];
}
Blas::Dgemv('N', nInteriorDofs, nInteriorDofs,
matsys->Scale(),
&((matsys->GetOwnedMatrix())->GetPtr())[0],
nInteriorDofs,rhs.get(),1,0.0,result.get(),1);
for (i = 0; i < nInteriorDofs; i++)
{
outarray[ mapArray[i] ] = result[i];
}
}
}
}
void Nektar::LocalRegions::QuadExp::v_GeneralMatrixOp_MatOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdRegions::StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdQuadExp.

Definition at line 2064 of file QuadExp.cpp.

References Nektar::LocalRegions::Expansion::GetLocMatrix(), Nektar::StdRegions::StdExpansion::m_ncoeffs, and Vmath::Vcopy().

{
MatrixKey newkey(mkey);
if (inarray.get() == outarray.get())
{
Array<OneD,NekDouble> tmp(m_ncoeffs);
Vmath::Vcopy(m_ncoeffs,inarray.get(),1,tmp.get(),1);
Blas::Dgemv('N',m_ncoeffs,m_ncoeffs, mat->Scale(),
(mat->GetOwnedMatrix())->GetPtr().get(), m_ncoeffs,
tmp.get(), 1, 0.0, outarray.get(), 1);
}
else
{
Blas::Dgemv('N',m_ncoeffs,m_ncoeffs,mat->Scale(),
(mat->GetOwnedMatrix())->GetPtr().get(), m_ncoeffs,
inarray.get(), 1, 0.0, outarray.get(), 1);
}
}
DNekMatSharedPtr Nektar::LocalRegions::QuadExp::v_GenMatrix ( const StdRegions::StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdQuadExp.

Definition at line 1501 of file QuadExp.cpp.

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

const LibUtilities::BasisSharedPtr & Nektar::LocalRegions::QuadExp::v_GetBasis ( int  dir) const
protectedvirtual

Definition at line 1489 of file QuadExp.cpp.

References ASSERTL1, and Nektar::StdRegions::StdExpansion::m_base.

{
ASSERTL1(dir >= 0 &&dir <= 1, "input dir is out of range");
return m_base[dir];
}
StdRegions::Orientation Nektar::LocalRegions::QuadExp::v_GetCartesianEorient ( int  edge)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1483 of file QuadExp.cpp.

References Nektar::LocalRegions::Expansion2D::GetGeom2D().

{
return GetGeom2D()->GetCartesianEorient(edge);
}
void Nektar::LocalRegions::QuadExp::v_GetCoord ( const Array< OneD, const NekDouble > &  Lcoords,
Array< OneD, NekDouble > &  coords 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 604 of file QuadExp.cpp.

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

{
int i;
ASSERTL1(Lcoords[0] >= -1.0 && Lcoords[1] <= 1.0 &&
Lcoords[1] >= -1.0 && Lcoords[1] <=1.0,
"Local coordinates are not in region [-1,1]");
m_geom->FillGeom();
for (i = 0; i < m_geom->GetCoordim(); ++i)
{
coords[i] = m_geom->GetCoord(i,Lcoords);
}
}
int Nektar::LocalRegions::QuadExp::v_GetCoordim ( void  )
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion2D.

Definition at line 1404 of file QuadExp.cpp.

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

{
return m_geom->GetCoordim();
}
void Nektar::LocalRegions::QuadExp::v_GetCoords ( Array< OneD, NekDouble > &  coords_1,
Array< OneD, NekDouble > &  coords_2,
Array< OneD, NekDouble > &  coords_3 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdQuadExp.

Definition at line 595 of file QuadExp.cpp.

{
Expansion::v_GetCoords(coords_0, coords_1, coords_2);
}
void Nektar::LocalRegions::QuadExp::v_GetEdgeInterpVals ( const int  edge,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 791 of file QuadExp.cpp.

References ASSERTL0, Vmath::Ddot(), Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::StdRegions::eFactorGaussEdge, Nektar::StdRegions::eInterpGauss, Nektar::StdRegions::StdExpansion::m_base, and m_matrixManager.

Referenced by v_ComputeEdgeNormal(), v_GetEdgePhysVals(), and v_GetEdgeQFactors().

{
int i;
int nq0 = m_base[0]->GetNumPoints();
int nq1 = m_base[1]->GetNumPoints();
StdRegions::StdMatrixKey key(
DetShapeType(),*this,factors);
switch (edge)
{
case 0:
{
for (i = 0; i < nq0; i++)
{
outarray[i] = Blas::Ddot(
nq1, mat_gauss->GetOwnedMatrix()->GetPtr().get(),
1, &inarray[i], nq0);
}
break;
}
case 1:
{
for (i = 0; i < nq1; i++)
{
outarray[i] = Blas::Ddot(
nq0, mat_gauss->GetOwnedMatrix()->GetPtr().get(),
1, &inarray[i * nq0], 1);
}
break;
}
case 2:
{
for (i = 0; i < nq0; i++)
{
outarray[i] = Blas::Ddot(
nq1, mat_gauss->GetOwnedMatrix()->GetPtr().get(),
1, &inarray[i], nq0);
}
break;
}
case 3:
{
for (i = 0; i < nq1; i++)
{
outarray[i] = Blas::Ddot(
nq0, mat_gauss->GetOwnedMatrix()->GetPtr().get(),
1, &inarray[i * nq0], 1);
}
break;
}
default:
ASSERTL0(false,"edge value (< 3) is out of range");
break;
}
}
void Nektar::LocalRegions::QuadExp::v_GetEdgePhysVals ( const int  edge,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Extract the physical values along edge edge from inarray into outarray following the local edge orientation and point distribution defined by defined in EdgeExp.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 651 of file QuadExp.cpp.

References ASSERTL0, Nektar::StdRegions::eForwards, Nektar::StdRegions::StdExpansion::GetEorient(), Nektar::StdRegions::StdExpansion::m_base, and Vmath::Vcopy().

Referenced by v_GetTracePhysVals().

{
int nquad0 = m_base[0]->GetNumPoints();
int nquad1 = m_base[1]->GetNumPoints();
switch(edge)
{
case 0:
if (edgedir == StdRegions::eForwards)
{
Vmath::Vcopy(nquad0,&(inarray[0]),1,&(outarray[0]),1);
}
else
{
Vmath::Vcopy(nquad0,&(inarray[0])+(nquad0-1),-1,
&(outarray[0]),1);
}
break;
case 1:
if (edgedir == StdRegions::eForwards)
{
Vmath::Vcopy(nquad1,&(inarray[0])+(nquad0-1),nquad0,
&(outarray[0]),1);
}
else
{
Vmath::Vcopy(nquad1,&(inarray[0])+(nquad0*nquad1-1),
-nquad0, &(outarray[0]),1);
}
break;
case 2:
if (edgedir == StdRegions::eForwards)
{
Vmath::Vcopy(nquad0,&(inarray[0])+(nquad0*nquad1-1),-1,
&(outarray[0]),1);
}
else
{
Vmath::Vcopy(nquad0,&(inarray[0])+nquad0*(nquad1-1),1,
&(outarray[0]),1);
}
break;
case 3:
if (edgedir == StdRegions::eForwards)
{
Vmath::Vcopy(nquad1,&(inarray[0]) + nquad0*(nquad1-1),
-nquad0,&(outarray[0]),1);
}
else
{
Vmath::Vcopy(nquad1,&(inarray[0]),nquad0,
&(outarray[0]),1);
}
break;
default:
ASSERTL0(false,"edge value (< 3) is out of range");
break;
}
}
void Nektar::LocalRegions::QuadExp::v_GetEdgePhysVals ( const int  edge,
const StdRegions::StdExpansionSharedPtr EdgeExp,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Definition at line 727 of file QuadExp.cpp.

References ASSERTL0, Nektar::StdRegions::eBackwards, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::StdRegions::StdExpansion::GetCartesianEorient(), Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::LibUtilities::Interp1D(), Nektar::StdRegions::StdExpansion::m_base, Vmath::Reverse(), v_GetEdgeInterpVals(), and Vmath::Vcopy().

{
int nquad0 = m_base[0]->GetNumPoints();
int nquad1 = m_base[1]->GetNumPoints();
// Implementation for all the basis except Gauss points
if (m_base[0]->GetPointsType() !=
{
// get points in Cartesian orientation
switch (edge)
{
case 0:
Vmath::Vcopy(nquad0,&(inarray[0]),1,&(outarray[0]),1);
break;
case 1:
Vmath::Vcopy(nquad1,&(inarray[0])+(nquad0-1),
nquad0,&(outarray[0]),1);
break;
case 2:
Vmath::Vcopy(nquad0,&(inarray[0])+nquad0*(nquad1-1),1,
&(outarray[0]),1);
break;
case 3:
Vmath::Vcopy(nquad1,&(inarray[0]),nquad0,
&(outarray[0]),1);
break;
default:
ASSERTL0(false,"edge value (< 3) is out of range");
break;
}
}
else
{
QuadExp::v_GetEdgeInterpVals(edge, inarray, outarray);
}
// Interpolate if required
if (m_base[edge%2]->GetPointsKey() !=
EdgeExp->GetBasis(0)->GetPointsKey())
{
Array<OneD,NekDouble> outtmp(max(nquad0,nquad1));
outtmp = outarray;
m_base[edge%2]->GetPointsKey(),outtmp,
EdgeExp->GetBasis(0)->GetPointsKey(),outarray);
}
//Reverse data if necessary
{
Vmath::Reverse(EdgeExp->GetNumPoints(0),&outarray[0],1,
&outarray[0],1);
}
}
void Nektar::LocalRegions::QuadExp::v_GetEdgeQFactors ( const int  edge,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 856 of file QuadExp.cpp.

References ASSERTL0, Nektar::SpatialDomains::eDeformed, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_metricinfo, Vmath::Reverse(), v_GetEdgeInterpVals(), Vmath::Vcopy(), and Vmath::Vmul().

{
int i;
int nquad0 = m_base[0]->GetNumPoints();
int nquad1 = m_base[1]->GetNumPoints();
const Array<OneD, const NekDouble>& jac = m_metricinfo->GetJac(ptsKeys);
const Array<TwoD, const NekDouble>& df = m_metricinfo->GetDerivFactors(ptsKeys);
Array<OneD, NekDouble> j (max(nquad0, nquad1), 0.0);
Array<OneD, NekDouble> g0(max(nquad0, nquad1), 0.0);
Array<OneD, NekDouble> g1(max(nquad0, nquad1), 0.0);
Array<OneD, NekDouble> g2(max(nquad0, nquad1), 0.0);
Array<OneD, NekDouble> g3(max(nquad0, nquad1), 0.0);
{
// Implementation for all the basis except Gauss points
&& m_base[1]->GetPointsType() !=
{
switch (edge)
{
case 0:
Vmath::Vcopy(nquad0, &(df[1][0]),
1, &(g1[0]), 1);
Vmath::Vcopy(nquad0, &(df[3][0]),
1, &(g3[0]), 1);
Vmath::Vcopy(nquad0, &(jac[0]),1, &(j[0]), 1);
for (i = 0; i < nquad0; ++i)
{
outarray[i] = j[i]*sqrt(g1[i]*g1[i]
+ g3[i]*g3[i]);
}
break;
case 1:
Vmath::Vcopy(nquad1,
&(df[0][0])+(nquad0-1), nquad0,
&(g0[0]), 1);
Vmath::Vcopy(nquad1,
&(df[2][0])+(nquad0-1), nquad0,
&(g2[0]), 1);
Vmath::Vcopy(nquad1,
&(jac[0])+(nquad0-1), nquad0,
&(j[0]), 1);
for (i = 0; i < nquad0; ++i)
{
outarray[i] = j[i]*sqrt(g0[i]*g0[i] +
g2[i]*g2[i]);
}
break;
case 2:
Vmath::Vcopy(nquad0,
&(df[1][0])+(nquad0*nquad1-1), -1,
&(g1[0]), 1);
Vmath::Vcopy(nquad0,
&(df[3][0])+(nquad0*nquad1-1), -1,
&(g3[0]), 1);
Vmath::Vcopy(nquad0,
&(jac[0])+(nquad0*nquad1-1), -1,
&(j[0]), 1);
for (i = 0; i < nquad0; ++i)
{
outarray[i] = j[i]*sqrt(g1[i]*g1[i]
+ g3[i]*g3[i]);
}
break;
case 3:
Vmath::Vcopy(nquad1,
&(df[0][0])+nquad0*(nquad1-1),
-nquad0,&(g0[0]), 1);
Vmath::Vcopy(nquad1,
&(df[2][0])+nquad0*(nquad1-1),
-nquad0,&(g2[0]), 1);
Vmath::Vcopy(nquad1,
&(jac[0])+nquad0*(nquad1-1), -nquad0,
&(j[0]), 1);
for (i = 0; i < nquad0; ++i)
{
outarray[i] = j[i]*sqrt(g0[i]*g0[i] +
g2[i]*g2[i]);
}
break;
default:
ASSERTL0(false,"edge value (< 3) is out of range");
break;
}
}
else
{
int nqtot = nquad0 * nquad1;
Array<OneD, NekDouble> tmp_gmat0(nqtot, 0.0);
Array<OneD, NekDouble> tmp_gmat1(nqtot, 0.0);
Array<OneD, NekDouble> tmp_gmat2(nqtot, 0.0);
Array<OneD, NekDouble> tmp_gmat3(nqtot, 0.0);
Array<OneD, NekDouble> g0_edge(max(nquad0, nquad1), 0.0);
Array<OneD, NekDouble> g1_edge(max(nquad0, nquad1), 0.0);
Array<OneD, NekDouble> g2_edge(max(nquad0, nquad1), 0.0);
Array<OneD, NekDouble> g3_edge(max(nquad0, nquad1), 0.0);
Array<OneD, NekDouble> jac_edge(max(nquad0, nquad1), 0.0);
switch (edge)
{
case 0:
Vmath::Vmul(nqtot,&(df[1][0]),1,&jac[0],1,
&(tmp_gmat1[0]),1);
Vmath::Vmul(nqtot,&(df[3][0]),1,&jac[0],1,
&(tmp_gmat3[0]),1);
edge, tmp_gmat1, g1_edge);
edge, tmp_gmat3, g3_edge);
for (i = 0; i < nquad0; ++i)
{
outarray[i] = sqrt(g1_edge[i]*g1_edge[i] +
g3_edge[i]*g3_edge[i]);
}
break;
case 1:
Vmath::Vmul(nqtot,
&(df[0][0]), 1,
&jac[0], 1,
&(tmp_gmat0[0]), 1);
Vmath::Vmul(nqtot,
&(df[2][0]), 1,
&jac[0], 1,
&(tmp_gmat2[0]),
1);
edge, tmp_gmat0, g0_edge);
edge, tmp_gmat2, g2_edge);
for (i = 0; i < nquad1; ++i)
{
outarray[i] = sqrt(g0_edge[i]*g0_edge[i]
+ g2_edge[i]*g2_edge[i]);
}
break;
case 2:
Vmath::Vmul(nqtot,
&(df[1][0]), 1,
&jac[0], 1,
&(tmp_gmat1[0]), 1);
Vmath::Vmul(nqtot,
&(df[3][0]), 1,
&jac[0], 1,
&(tmp_gmat3[0]),1);
edge, tmp_gmat1, g1_edge);
edge, tmp_gmat3, g3_edge);
for (i = 0; i < nquad0; ++i)
{
outarray[i] = sqrt(g1_edge[i]*g1_edge[i]
+ g3_edge[i]*g3_edge[i]);
}
Vmath::Reverse(nquad0,&outarray[0],1,&outarray[0],1);
break;
case 3:
Vmath::Vmul(nqtot,
&(df[0][0]), 1,
&jac[0], 1,
&(tmp_gmat0[0]), 1);
Vmath::Vmul(nqtot,
&(df[2][0]),1,
&jac[0], 1,
&(tmp_gmat2[0]),1);
edge, tmp_gmat0, g0_edge);
edge, tmp_gmat2, g2_edge);
for (i = 0; i < nquad1; ++i)
{
outarray[i] = sqrt(g0_edge[i]*g0_edge[i] +
g2_edge[i]*g2_edge[i]);
}
&outarray[0], 1,
&outarray[0], 1);
break;
default:
ASSERTL0(false,"edge value (< 3) is out of range");
break;
}
}
}
else
{
switch (edge)
{
case 0:
for (i = 0; i < nquad0; ++i)
{
outarray[i] = jac[0]*sqrt(df[1][0]*df[1][0] +
df[3][0]*df[3][0]);
}
break;
case 1:
for (i = 0; i < nquad1; ++i)
{
outarray[i] = jac[0]*sqrt(df[0][0]*df[0][0] +
df[2][0]*df[2][0]);
}
break;
case 2:
for (i = 0; i < nquad0; ++i)
{
outarray[i] = jac[0]*sqrt(df[1][0]*df[1][0] +
df[3][0]*df[3][0]);
}
break;
case 3:
for (i = 0; i < nquad1; ++i)
{
outarray[i] = jac[0]*sqrt(df[0][0]*df[0][0] +
df[2][0]*df[2][0]);
}
break;
default:
ASSERTL0(false,"edge value (< 3) is out of range");
break;
}
}
}
StdRegions::Orientation Nektar::LocalRegions::QuadExp::v_GetEorient ( int  edge)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1477 of file QuadExp.cpp.

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

{
return m_geom->GetEorient(edge);
}
DNekScalMatSharedPtr Nektar::LocalRegions::QuadExp::v_GetLocMatrix ( const MatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1977 of file QuadExp.cpp.

References m_matrixManager.

{
return m_matrixManager[mkey];
}
DNekScalBlkMatSharedPtr Nektar::LocalRegions::QuadExp::v_GetLocStaticCondMatrix ( const MatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1983 of file QuadExp.cpp.

References m_staticCondMatrixManager.

{
}
const SpatialDomains::GeomFactorsSharedPtr & Nektar::LocalRegions::QuadExp::v_GetMetricInfo ( ) const
protectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1398 of file QuadExp.cpp.

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

{
return m_metricinfo;
}
int Nektar::LocalRegions::QuadExp::v_GetNumPoints ( const int  dir) const
protectedvirtual

Definition at line 1496 of file QuadExp.cpp.

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

{
return GetNumPoints(dir);
}
void Nektar::LocalRegions::QuadExp::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 716 of file QuadExp.cpp.

References v_GetEdgePhysVals().

{
v_GetEdgePhysVals(edge,EdgeExp,inarray,outarray);
}
void Nektar::LocalRegions::QuadExp::v_HelmholtzMatrixOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdRegions::StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdQuadExp.

Definition at line 2055 of file QuadExp.cpp.

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

{
QuadExp::HelmholtzMatrixOp_MatFree(inarray, outarray, mkey);
}
NekDouble Nektar::LocalRegions::QuadExp::v_Integral ( const Array< OneD, const NekDouble > &  inarray)
protectedvirtual

Integrates the specified function over the domain.

See Also
StdRegions::StdExpansion::Integral.

Reimplemented from Nektar::StdRegions::StdQuadExp.

Definition at line 83 of file QuadExp.cpp.

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

{
int nquad0 = m_base[0]->GetNumPoints();
int nquad1 = m_base[1]->GetNumPoints();
Array<OneD, const NekDouble> jac = m_metricinfo->GetJac(GetPointsKeys());
NekDouble ival;
Array<OneD,NekDouble> tmp(nquad0*nquad1);
// multiply inarray with Jacobian
{
Vmath::Vmul(nquad0*nquad1, jac, 1, inarray, 1, tmp, 1);
}
else
{
Vmath::Smul(nquad0*nquad1, jac[0], inarray, 1, tmp, 1);
}
// call StdQuadExp version;
return ival;
}
void Nektar::LocalRegions::QuadExp::v_IProductWRTBase ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Calculate the inner product of inarray with respect to the basis B=base0*base1 and put into outarray.

$ \begin{array}{rcl} I_{pq} = (\phi_q \phi_q, u) & = & \sum_{i=0}^{nq_0} \sum_{j=0}^{nq_1} \phi_p(\xi_{0,i}) \phi_q(\xi_{1,j}) w^0_i w^1_j u(\xi_{0,i} \xi_{1,j}) \\ & = & \sum_{i=0}^{nq_0} \phi_p(\xi_{0,i}) \sum_{j=0}^{nq_1} \phi_q(\xi_{1,j}) \tilde{u}_{i,j} \end{array} $

where

$ \tilde{u}_{i,j} = w^0_i w^1_j u(\xi_{0,i},\xi_{1,j}) $

which can be implemented as

$ f_{qi} = \sum_{j=0}^{nq_1} \phi_q(\xi_{1,j}) \tilde{u}_{i,j} = {\bf B_1 U} $ $ I_{pq} = \sum_{i=0}^{nq_0} \phi_p(\xi_{0,i}) f_{qi} = {\bf B_0 F} $

Reimplemented from Nektar::StdRegions::StdQuadExp.

Definition at line 409 of file QuadExp.cpp.

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

{
if (m_base[0]->Collocation() && m_base[1]->Collocation())
{
MultiplyByQuadratureMetric(inarray,outarray);
}
else
{
IProductWRTBase_SumFac(inarray,outarray);
}
}
void Nektar::LocalRegions::QuadExp::v_IProductWRTBase_MatOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdQuadExp.

Definition at line 452 of file QuadExp.cpp.

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

{
int nq = GetTotPoints();
MatrixKey
DNekScalMatSharedPtr iprodmat = m_matrixManager[iprodmatkey];
Blas::Dgemv('N',m_ncoeffs,nq,iprodmat->Scale(),
(iprodmat->GetOwnedMatrix())->GetPtr().get(),
m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
}
void Nektar::LocalRegions::QuadExp::v_IProductWRTBase_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdQuadExp.

Definition at line 433 of file QuadExp.cpp.

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

{
int nquad0 = m_base[0]->GetNumPoints();
int nquad1 = m_base[1]->GetNumPoints();
int order0 = m_base[0]->GetNumModes();
Array<OneD,NekDouble> tmp(nquad0*nquad1+nquad1*order0);
Array<OneD,NekDouble> wsp(tmp+nquad0*nquad1);
m_base[0]->GetBdata(),
m_base[1]->GetBdata(),
tmp,outarray,wsp,true,true);
}
void Nektar::LocalRegions::QuadExp::v_IProductWRTDerivBase ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdQuadExp.

Definition at line 424 of file QuadExp.cpp.

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

{
IProductWRTDerivBase_SumFac(dir,inarray,outarray);
}
void Nektar::LocalRegions::QuadExp::v_IProductWRTDerivBase_MatOp ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdQuadExp.

Definition at line 524 of file QuadExp.cpp.

References ASSERTL1, Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::StdRegions::eIProductWRTDerivBase0, Nektar::StdRegions::eIProductWRTDerivBase1, Nektar::StdRegions::eIProductWRTDerivBase2, Nektar::StdRegions::StdExpansion::GetTotPoints(), m_matrixManager, and Nektar::StdRegions::StdExpansion::m_ncoeffs.

{
int nq = GetTotPoints();
switch (dir)
{
case 0:
{
}
break;
case 1:
{
}
break;
case 2:
{
}
break;
default:
{
ASSERTL1(false,"input dir is out of range");
}
break;
}
MatrixKey iprodmatkey(mtype,DetShapeType(),*this);
DNekScalMatSharedPtr iprodmat = m_matrixManager[iprodmatkey];
Blas::Dgemv('N', m_ncoeffs, nq, iprodmat->Scale(),
(iprodmat->GetOwnedMatrix())->GetPtr().get(),
m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
}
void Nektar::LocalRegions::QuadExp::v_IProductWRTDerivBase_SumFac ( const int  dir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdQuadExp.

Definition at line 468 of file QuadExp.cpp.

References ASSERTL1, Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion2D::IProductWRTBase_SumFacKernel(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_geom, Nektar::LocalRegions::Expansion::m_metricinfo, Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::StdRegions::StdExpansion::MultiplyByQuadratureMetric(), Vmath::Smul(), Vmath::Vadd(), and Vmath::Vmul().

{
ASSERTL1((dir==0) || (dir==1) || (dir==2),
"Invalid direction.");
ASSERTL1((dir==2) ? (m_geom->GetCoordim() ==3):true,
"Invalid direction.");
int nquad0 = m_base[0]->GetNumPoints();
int nquad1 = m_base[1]->GetNumPoints();
int nqtot = nquad0*nquad1;
int nmodes0 = m_base[0]->GetNumModes();
const Array<TwoD, const NekDouble>& df = m_metricinfo->GetDerivFactors(GetPointsKeys());
Array<OneD, NekDouble> tmp1(2*nqtot+m_ncoeffs+nmodes0*nquad1);
Array<OneD, NekDouble> tmp2(tmp1 + nqtot);
Array<OneD, NekDouble> tmp3(tmp1 + 2*nqtot);
Array<OneD, NekDouble> tmp4(tmp1 + 2*nqtot+m_ncoeffs);
{
Vmath::Vmul(nqtot,
&df[2*dir][0], 1,
inarray.get(), 1,
tmp1.get(), 1);
Vmath::Vmul(nqtot,
&df[2*dir+1][0], 1,
inarray.get(), 1,
tmp2.get(),1);
}
else
{
Vmath::Smul(nqtot,
df[2*dir][0], inarray.get(), 1,
tmp1.get(), 1);
Vmath::Smul(nqtot,
df[2*dir+1][0], inarray.get(), 1,
tmp2.get(), 1);
}
m_base[0]->GetDbdata(), m_base[1]->GetBdata(),
tmp1, tmp3, tmp4, false, true);
m_base[0]->GetBdata() , m_base[1]->GetDbdata(),
tmp2, outarray, tmp4, true, false);
Vmath::Vadd(m_ncoeffs, tmp3, 1, outarray, 1, outarray, 1);
}
void Nektar::LocalRegions::QuadExp::v_LaplacianMatrixOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdRegions::StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdQuadExp.

Definition at line 2004 of file QuadExp.cpp.

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

{
QuadExp::LaplacianMatrixOp_MatFree(inarray, outarray, mkey);
}
void Nektar::LocalRegions::QuadExp::v_LaplacianMatrixOp ( const int  k1,
const int  k2,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdRegions::StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdQuadExp.

Definition at line 2013 of file QuadExp.cpp.

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

{
k1, k2, inarray, outarray, mkey);
}
void Nektar::LocalRegions::QuadExp::v_LaplacianMatrixOp_MatFree_Kernel ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
Array< OneD, NekDouble > &  wsp 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2139 of file QuadExp.cpp.

References ASSERTL1, Nektar::LocalRegions::Expansion::ComputeLaplacianMetric(), Nektar::StdRegions::StdExpansion2D::IProductWRTBase_SumFacKernel(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_metrics, Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::LocalRegions::MetricLaplacian00, Nektar::LocalRegions::MetricLaplacian01, Nektar::LocalRegions::MetricLaplacian11, Nektar::StdRegions::StdExpansion2D::PhysTensorDeriv(), Vmath::Vadd(), and Vmath::Vvtvvtp().

{
if (m_metrics.count(MetricLaplacian00) == 0)
{
}
int nquad0 = m_base[0]->GetNumPoints();
int nquad1 = m_base[1]->GetNumPoints();
int nqtot = nquad0*nquad1;
int nmodes0 = m_base[0]->GetNumModes();
int nmodes1 = m_base[1]->GetNumModes();
int wspsize = max(max(max(nqtot,m_ncoeffs),nquad1*nmodes0),nquad0*nmodes1);
ASSERTL1(wsp.num_elements() >= 3*wspsize,
"Workspace is of insufficient size.");
const Array<OneD, const NekDouble>& base0 = m_base[0]->GetBdata();
const Array<OneD, const NekDouble>& base1 = m_base[1]->GetBdata();
const Array<OneD, const NekDouble>& dbase0 = m_base[0]->GetDbdata();
const Array<OneD, const NekDouble>& dbase1 = m_base[1]->GetDbdata();
const Array<OneD, const NekDouble>& metric00 = m_metrics[MetricLaplacian00];
const Array<OneD, const NekDouble>& metric01 = m_metrics[MetricLaplacian01];
const Array<OneD, const NekDouble>& metric11 = m_metrics[MetricLaplacian11];
// Allocate temporary storage
Array<OneD,NekDouble> wsp0(wsp);
Array<OneD,NekDouble> wsp1(wsp+wspsize);
Array<OneD,NekDouble> wsp2(wsp+2*wspsize);
StdExpansion2D::PhysTensorDeriv(inarray,wsp1,wsp2);
// wsp0 = k = g0 * wsp1 + g1 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
// wsp2 = l = g1 * wsp1 + g2 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
// where g0, g1 and g2 are the metric terms set up in the GeomFactors class
// especially for this purpose
Vmath::Vvtvvtp(nqtot,&metric00[0],1,&wsp1[0],1,&metric01[0],1,&wsp2[0],1,&wsp0[0],1);
Vmath::Vvtvvtp(nqtot,&metric01[0],1,&wsp1[0],1,&metric11[0],1,&wsp2[0],1,&wsp2[0],1);
// outarray = m = (D_xi1 * B)^T * k
// wsp1 = n = (D_xi2 * B)^T * l
IProductWRTBase_SumFacKernel(dbase0,base1,wsp0,outarray,wsp1,false,true);
IProductWRTBase_SumFacKernel(base0,dbase1,wsp2,wsp1, wsp0,true,false);
// outarray = outarray + wsp1
// = L * u_hat
Vmath::Vadd(m_ncoeffs,wsp1.get(),1,outarray.get(),1,outarray.get(),1);
}
void Nektar::LocalRegions::QuadExp::v_MassLevelCurvatureMatrixOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdRegions::StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2045 of file QuadExp.cpp.

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

{
inarray, outarray, mkey);
}
void Nektar::LocalRegions::QuadExp::v_MassMatrixOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdRegions::StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdQuadExp.

Definition at line 1995 of file QuadExp.cpp.

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

{
StdExpansion::MassMatrixOp_MatFree(inarray, outarray, mkey);
}
void Nektar::LocalRegions::QuadExp::v_NormVectorIProductWRTBase ( const Array< OneD, const NekDouble > &  Fx,
const Array< OneD, const NekDouble > &  Fy,
const Array< OneD, const NekDouble > &  Fz,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 565 of file QuadExp.cpp.

References Nektar::SpatialDomains::eDeformed, Nektar::LocalRegions::Expansion2D::GetLeftAdjacentElementExp(), Nektar::LocalRegions::Expansion2D::GetLeftAdjacentElementFace(), Nektar::StdRegions::StdExpansion::IProductWRTBase(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_metricinfo, Vmath::Svtsvtp(), Vmath::Svtvp(), Vmath::Vvtvp(), and Vmath::Vvtvvtp().

{
int nq = m_base[0]->GetNumPoints()*m_base[1]->GetNumPoints();
Array<OneD, NekDouble> Fn(nq);
const Array<OneD, const Array<OneD, NekDouble> > &normals =
GetLeftAdjacentElementExp()->GetFaceNormal(
{
Vmath::Vvtvvtp(nq,&normals[0][0],1,&Fx[0],1,
&normals[1][0],1,&Fy[0],1,&Fn[0],1);
Vmath::Vvtvp (nq,&normals[2][0],1,&Fz[0],1,&Fn[0],1,&Fn[0],1);
}
else
{
Vmath::Svtsvtp(nq,normals[0][0],&Fx[0],1,
normals[1][0],&Fy[0],1,&Fn[0],1);
Vmath::Svtvp (nq,normals[2][0],&Fz[0],1,&Fn[0],1,&Fn[0],1);
}
IProductWRTBase(Fn,outarray);
}
void Nektar::LocalRegions::QuadExp::v_PhysDeriv ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d0,
Array< OneD, NekDouble > &  out_d1,
Array< OneD, NekDouble > &  out_d2 = NullNekDouble1DArray 
)
protectedvirtual

Calculate the derivative of the physical points.

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

Reimplemented from Nektar::StdRegions::StdQuadExp.

Definition at line 110 of file QuadExp.cpp.

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

Referenced by v_PhysDeriv(), and v_PhysDirectionalDeriv().

{
int nquad0 = m_base[0]->GetNumPoints();
int nquad1 = m_base[1]->GetNumPoints();
int nqtot = nquad0*nquad1;
const Array<TwoD, const NekDouble>& df = m_metricinfo->GetDerivFactors(GetPointsKeys());
Array<OneD,NekDouble> diff0(2*nqtot);
Array<OneD,NekDouble> diff1(diff0+nqtot);
StdQuadExp::v_PhysDeriv(inarray, diff0, diff1);
{
if (out_d0.num_elements())
{
Vmath::Vmul (nqtot, df[0], 1, diff0, 1, out_d0, 1);
Vmath::Vvtvp (nqtot, df[1], 1, diff1, 1, out_d0, 1,
out_d0,1);
}
if(out_d1.num_elements())
{
Vmath::Vmul (nqtot,df[2],1,diff0,1, out_d1, 1);
Vmath::Vvtvp (nqtot,df[3],1,diff1,1, out_d1, 1, out_d1,1);
}
if (out_d2.num_elements())
{
Vmath::Vmul (nqtot,df[4],1,diff0,1, out_d2, 1);
Vmath::Vvtvp (nqtot,df[5],1,diff1,1, out_d2, 1, out_d2,1);
}
}
else // regular geometry
{
if (out_d0.num_elements())
{
Vmath::Smul (nqtot, df[0][0], diff0, 1, out_d0, 1);
Blas::Daxpy (nqtot, df[1][0], diff1, 1, out_d0, 1);
}
if (out_d1.num_elements())
{
Vmath::Smul (nqtot, df[2][0], diff0, 1, out_d1, 1);
Blas::Daxpy (nqtot, df[3][0], diff1, 1, out_d1, 1);
}
if (out_d2.num_elements())
{
Vmath::Smul (nqtot, df[4][0], diff0, 1, out_d2, 1);
Blas::Daxpy (nqtot, df[5][0], diff1, 1, out_d2, 1);
}
}
}
void Nektar::LocalRegions::QuadExp::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::StdQuadExp.

Definition at line 169 of file QuadExp.cpp.

References ASSERTL1, Nektar::NullNekDouble1DArray, and v_PhysDeriv().

{
switch (dir)
{
case 0:
{
v_PhysDeriv(inarray, outarray, NullNekDouble1DArray,
}
break;
case 1:
{
v_PhysDeriv(inarray, NullNekDouble1DArray, outarray,
}
break;
case 2:
{
}
break;
default:
{
ASSERTL1(false,"input dir is out of range");
}
break;
}
}
void Nektar::LocalRegions::QuadExp::v_PhysDirectionalDeriv ( const Array< OneD, const NekDouble > &  inarray,
const Array< OneD, const NekDouble > &  direction,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Physical derivative along a direction vector.

See Also
StdRegions::StdExpansion::PhysDirectionalDeriv

D_v = d/dx_v^s + d/dx_v^r

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 203 of file QuadExp.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, v_PhysDeriv(), Vmath::Vmul(), and Vmath::Vvtvp().

{
int nquad0 = m_base[0]->GetNumPoints();
int nquad1 = m_base[1]->GetNumPoints();
int nqtot = nquad0*nquad1;
const Array<TwoD, const NekDouble>& df = m_metricinfo->GetDerivFactors(GetPointsKeys());
Array<OneD,NekDouble> diff0(2*nqtot);
Array<OneD,NekDouble> diff1(diff0+nqtot);
StdQuadExp::v_PhysDeriv(inarray, diff0, diff1);
{
Array<OneD, Array<OneD, NekDouble> > tangmat(2);
// d/dx_v^s = v_x*ds/dx + v_y*ds/dy + v_z*dx/dz
for (int i=0; i< 2; ++i)
{
tangmat[i] = Array<OneD, NekDouble>(nqtot,0.0);
for (int k=0; k<(m_geom->GetCoordim()); ++k)
{
Vmath::Vvtvp(nqtot,
&df[2*k+i][0], 1,
&direction[k*nqtot], 1,
&tangmat[i][0], 1,
&tangmat[i][0], 1);
}
}
/// D_v = d/dx_v^s + d/dx_v^r
if (out.num_elements())
{
Vmath::Vmul (nqtot,
&tangmat[0][0], 1,
&diff0[0], 1,
&out[0], 1);
Vmath::Vvtvp (nqtot,
&tangmat[1][0], 1,
&diff1[0], 1,
&out[0], 1,
&out[0], 1);
}
}
else
{
ASSERTL1(m_metricinfo->GetGtype() ==
SpatialDomains::eDeformed,"Wrong route");
}
}
NekDouble Nektar::LocalRegions::QuadExp::v_PhysEvaluate ( const Array< OneD, const NekDouble > &  coords,
const Array< OneD, const NekDouble > &  physvals 
)
protectedvirtual

This function evaluates the expansion at a single (arbitrary) point of the domain.

This function is a wrapper around the virtual function v_PhysEvaluate()

Based on the value of the expansion at the quadrature points, this function calculates the value of the expansion at an arbitrary single points (with coordinates $ \mathbf{x_c}$ given by the pointer coords). This operation, equivalent to

\[ u(\mathbf{x_c}) = \sum_p \phi_p(\mathbf{x_c}) \hat{u}_p \]

is evaluated using Lagrangian interpolants through the quadrature points:

\[ u(\mathbf{x_c}) = \sum_p h_p(\mathbf{x_c}) u_p\]

This function requires that the physical value array $\mathbf{u}$ (implemented as the attribute #m_phys) is set.

Parameters
coordsthe coordinates of the single point
Returns
returns the value of the expansion at the single point

Reimplemented from Nektar::StdRegions::StdExpansion2D.

Definition at line 635 of file QuadExp.cpp.

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

Referenced by v_StdPhysEvaluate().

{
Array<OneD,NekDouble> Lcoord = Array<OneD, NekDouble>(2);
ASSERTL0(m_geom,"m_geom not defined");
m_geom->GetLocCoords(coord,Lcoord);
return StdQuadExp::v_PhysEvaluate(Lcoord, physvals);
}
void Nektar::LocalRegions::QuadExp::v_ReduceOrderCoeffs ( int  numMin,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdQuadExp.

Definition at line 2089 of file QuadExp.cpp.

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

{
int n_coeffs = inarray.num_elements();
Array<OneD, NekDouble> coeff (n_coeffs);
Array<OneD, NekDouble> coeff_tmp(n_coeffs, 0.0);
Array<OneD, NekDouble> tmp, tmp2;
int nmodes0 = m_base[0]->GetNumModes();
int nmodes1 = m_base[1]->GetNumModes();
int numMax = nmodes0;
Vmath::Vcopy(n_coeffs,inarray,1,coeff_tmp,1);
const LibUtilities::PointsKey Pkey0(
const LibUtilities::PointsKey Pkey1(
LibUtilities::BasisKey b0(
m_base[0]->GetBasisType(), nmodes0, Pkey0);
LibUtilities::BasisKey b1(
m_base[1]->GetBasisType(), nmodes1, Pkey1);
LibUtilities::BasisKey bortho0(
LibUtilities::eOrtho_A, nmodes0, Pkey0);
LibUtilities::BasisKey bortho1(
LibUtilities::eOrtho_A, nmodes1, Pkey1);
b0, b1, coeff_tmp, bortho0, bortho1, coeff);
Vmath::Zero(n_coeffs, coeff_tmp, 1);
int cnt = 0;
for (int i = 0; i < numMin+1; ++i)
{
Vmath::Vcopy(numMin,
tmp = coeff+cnt,1,
tmp2 = coeff_tmp+cnt,1);
cnt = i*numMax;
}
bortho0, bortho1, coeff_tmp,
b0, b1, outarray);
}
NekDouble Nektar::LocalRegions::QuadExp::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 627 of file QuadExp.cpp.

References v_PhysEvaluate().

{
// Evaluate point in local (eta) coordinates.
return StdQuadExp::v_PhysEvaluate(Lcoord,physvals);
}
void Nektar::LocalRegions::QuadExp::v_WeakDerivMatrixOp ( const int  i,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdRegions::StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdQuadExp.

Definition at line 2025 of file QuadExp.cpp.

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

{
StdExpansion::WeakDerivMatrixOp_MatFree(i, inarray, outarray, mkey);
}
void Nektar::LocalRegions::QuadExp::v_WeakDirectionalDerivMatrixOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdRegions::StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 2035 of file QuadExp.cpp.

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

Member Data Documentation

LibUtilities::NekManager<MatrixKey, DNekScalMat, MatrixKey::opLess> Nektar::LocalRegions::QuadExp::m_matrixManager
private

Definition at line 268 of file QuadExp.h.

Referenced by CreateMatrix(), v_FwdTrans(), v_GetEdgeInterpVals(), v_GetLocMatrix(), v_IProductWRTBase_MatOp(), and v_IProductWRTDerivBase_MatOp().

LibUtilities::NekManager<MatrixKey, DNekScalBlkMat, MatrixKey::opLess> Nektar::LocalRegions::QuadExp::m_staticCondMatrixManager
private

Definition at line 269 of file QuadExp.h.

Referenced by v_DropLocStaticCondMatrix(), v_FwdTrans_BndConstrained(), and v_GetLocStaticCondMatrix().