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

#include <PrismExp.h>

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

Public Member Functions

 PrismExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc, const SpatialDomains::PrismGeomSharedPtr &geom)
 Constructor using BasisKey class for quadrature points and order definition.
 PrismExp (const PrismExp &T)
 ~PrismExp ()
- Public Member Functions inherited from Nektar::StdRegions::StdPrismExp
 StdPrismExp ()
 StdPrismExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
 StdPrismExp (const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc, NekDouble *coeffs, NekDouble *phys)
 StdPrismExp (const StdPrismExp &T)
 ~StdPrismExp ()
- Public Member Functions inherited from Nektar::StdRegions::StdExpansion3D
 StdExpansion3D ()
 StdExpansion3D (int numcoeffs, const LibUtilities::BasisKey &Ba, const LibUtilities::BasisKey &Bb, const LibUtilities::BasisKey &Bc)
 StdExpansion3D (const StdExpansion3D &T)
virtual ~StdExpansion3D ()
void PhysTensorDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d1, Array< OneD, NekDouble > &outarray_d2, Array< OneD, NekDouble > &outarray_d3)
 Calculate the 3D derivative in the local tensor/collapsed coordinate at the physical points.
void BwdTrans_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
void IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
- 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.
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 GetFaceOrient (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)
void AddEdgeNormBoundaryInt (const int edge, boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
void AddEdgeNormBoundaryInt (const int edge, boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
void AddEdgeNormBoundaryBiInt (const int edge, boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
void AddFaceNormBoundaryInt (const int face, boost::shared_ptr< StdExpansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
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)
void DGDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, boost::shared_ptr< StdExpansion > > &EdgeExp, Array< OneD, Array< OneD, NekDouble > > &coeffs, Array< OneD, NekDouble > &outarray)
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_SetUpPhysNormals (const int edge)
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
virtual void v_NormVectorIProductWRTBase (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, const Array< OneD, const NekDouble > &Fz, Array< OneD, NekDouble > &outarray)
virtual StdRegions::Orientation v_GetEorient (int edge)
virtual StdRegions::Orientation v_GetCartesianEorient (int edge)
virtual StdRegions::Orientation v_GetPorient (int point)
NekDouble Linf (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &sol=NullNekDouble1DArray)
 Function to evaluate the discrete $ L_\infty$ error $ |\epsilon|_\infty = \max |u - u_{exact}|$ where $ u_{exact}$ is given by the array sol.
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)
template<class T >
boost::shared_ptr< T > as ()
- Public Member Functions inherited from Nektar::LocalRegions::Expansion3D
 Expansion3D (SpatialDomains::Geometry3DSharedPtr pGeom)
virtual ~Expansion3D ()
void SetFaceExp (const int face, Expansion2DSharedPtr &f)
Expansion2DSharedPtr GetFaceExp (const int face)
void SetTraceToGeomOrientation (Array< OneD, NekDouble > &inout)
 Align trace orientation with the geometry orientation.
void SetFaceToGeomOrientation (const int face, Array< OneD, NekDouble > &inout)
 Align face orientation with the geometry orientation.
void AddHDGHelmholtzFaceTerms (const NekDouble tau, const int edge, Array< OneD, NekDouble > &facePhys, const StdRegions::VarCoeffMap &dirForcing, Array< OneD, NekDouble > &outarray)
void AddNormTraceInt (const int dir, Array< OneD, StdRegions::StdExpansionSharedPtr > &FaceExp, Array< OneD, Array< OneD, NekDouble > > &faceCoeffs, Array< OneD, NekDouble > &outarray)
void AddNormTraceInt (const int dir, Array< OneD, const NekDouble > &inarray, Array< OneD, StdRegions::StdExpansionSharedPtr > &FaceExp, Array< OneD, NekDouble > &outarray, const StdRegions::VarCoeffMap &varcoeffs)
void AddFaceBoundaryInt (const int face, StdRegions::StdExpansionSharedPtr &FaceExp, Array< OneD, NekDouble > &facePhys, Array< OneD, NekDouble > &outarray, const StdRegions::VarCoeffMap &varcoeffs=StdRegions::NullVarCoeffMap)
SpatialDomains::Geometry3DSharedPtr GetGeom3D () 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
virtual const
SpatialDomains::GeomFactorsSharedPtr
v_GetMetricInfo () const
DNekMatSharedPtr BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
DNekMatSharedPtr BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd)

Protected Member Functions

virtual NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray)
 Integrate the physical point list inarray over prismatic region and return the value.
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)
 Calculate the derivative of the physical points.
virtual void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Forward transform from physical quadrature space stored in inarray and evaluate the expansion coefficients and store in (this)->m_coeffs.
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*base2 and put into outarray:
virtual void v_IProductWRTBase_SumFac (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Calculates the inner product $ I_{pqr} = (u, \partial_{x_i} \phi_{pqr}) $.
void v_IProductWRTDerivBase_SumFac (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_GetCoord (const Array< OneD, const NekDouble > &Lcoords, Array< OneD, NekDouble > &coords)
 Get the coordinates #coords at the local coordinates #Lcoords.
virtual void v_GetCoords (Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2, Array< OneD, NekDouble > &coords_3)
virtual NekDouble v_StdPhysEvaluate (const Array< OneD, const NekDouble > &Lcoord, const Array< OneD, const NekDouble > &physvals)
virtual NekDouble v_PhysEvaluate (const Array< OneD, const NekDouble > &coord, const Array< OneD, const NekDouble > &physvals)
 This function evaluates the expansion at a single (arbitrary) point of the domain.
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_GetFaceOrient (int face)
virtual void v_GetFacePhysVals (const int face, const StdRegions::StdExpansionSharedPtr &FaceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient)
virtual void v_GetTracePhysVals (const int face, const StdRegions::StdExpansionSharedPtr &FaceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient)
 Returns the physical values at the quadrature points of a face.
void v_ComputeFaceNormal (const int face)
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_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_ReduceOrderCoeffs (int numMin, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual DNekMatSharedPtr v_GenMatrix (const StdRegions::StdMatrixKey &mkey)
virtual DNekMatSharedPtr v_CreateStdMatrix (const StdRegions::StdMatrixKey &mkey)
virtual DNekScalMatSharedPtr v_GetLocMatrix (const MatrixKey &mkey)
virtual DNekScalBlkMatSharedPtr v_GetLocStaticCondMatrix (const MatrixKey &mkey)
void v_DropLocStaticCondMatrix (const MatrixKey &mkey)
DNekScalMatSharedPtr CreateMatrix (const MatrixKey &mkey)
DNekScalBlkMatSharedPtr CreateStaticCondMatrix (const MatrixKey &mkey)
- Protected Member Functions inherited from Nektar::StdRegions::StdPrismExp
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_StdPhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
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 > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
virtual void v_IProductWRTBase_MatOp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_IProductWRTBase_SumFacKernel (const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &base2, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1, bool doCheckCollDir2)
virtual void v_IProductWRTDerivBase_MatOp (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_LocCoordToLocCollapsed (const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
virtual void v_FillMode (const int mode, Array< OneD, NekDouble > &outarray)
virtual int v_GetNverts () const
virtual int v_GetNedges () const
virtual int v_GetNfaces () const
virtual LibUtilities::ShapeType v_DetShapeType () const
 Return Shape of region, using ShapeType enum list; i.e. prism.
virtual int v_NumBndryCoeffs () const
virtual int v_NumDGBndryCoeffs () const
virtual int v_GetEdgeNcoeffs (const int i) const
virtual int v_GetTotalEdgeIntNcoeffs () const
virtual int v_GetFaceNcoeffs (const int i) const
virtual int v_GetFaceNumPoints (const int i) const
virtual LibUtilities::PointsKey v_GetFacePointsKey (const int i, const int j) const
virtual const
LibUtilities::BasisKey 
v_DetFaceBasisKey (const int i, const int k) const
virtual int v_GetFaceIntNcoeffs (const int i) const
virtual int v_GetTotalFaceIntNcoeffs () const
virtual int v_CalcNumberOfCoefficients (const std::vector< unsigned int > &nummodes, int &modes_offset)
virtual LibUtilities::BasisType v_GetEdgeBasisType (const int i) const
virtual bool v_IsBoundaryInteriorExpansion ()
virtual void v_GetFaceToElementMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray, int nummodesA=-1, int nummodesB=-1)
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)
virtual void v_GetFaceInteriorMap (const int fid, const Orientation faceOrient, Array< OneD, unsigned int > &maparray, Array< OneD, int > &signarray)
virtual void v_GetInteriorMap (Array< OneD, unsigned int > &outarray)
virtual void v_GetBoundaryMap (Array< OneD, unsigned int > &outarray)
virtual void v_MultiplyByStdQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void v_SVVLaplacianFilter (Array< OneD, NekDouble > &array, const StdMatrixKey &mkey)
- Protected Member Functions inherited from Nektar::StdRegions::StdExpansion3D
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)
virtual void v_NegateFaceNormal (const int face)
- 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)
virtual void v_AddEdgeNormBoundaryInt (const int edge, boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
virtual void v_AddEdgeNormBoundaryInt (const int edge, boost::shared_ptr< StdExpansion > &EdgeExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
virtual void v_AddFaceNormBoundaryInt (const int face, boost::shared_ptr< StdExpansion > &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
- Protected Member Functions inherited from Nektar::LocalRegions::Expansion3D
virtual void v_DGDeriv (const int dir, const Array< OneD, const NekDouble > &incoeffs, Array< OneD, StdRegions::StdExpansionSharedPtr > &FaceExp, Array< OneD, Array< OneD, NekDouble > > &faceCoeffs, Array< OneD, NekDouble > &out_d)
 Evaluate coefficients of weak deriviative in the direction dir given the input coefficicents incoeffs and the imposed boundary values in EdgeExp (which will have its phys space updated).
virtual void v_AddFaceNormBoundaryInt (const int face, StdRegions::StdExpansionSharedPtr &FaceExp, const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
virtual void v_AddRobinMassMatrix (const int face, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat)
virtual NekDouble v_Integrate (const Array< OneD, const NekDouble > &inarray)
virtual Array< OneD, unsigned int > v_GetEdgeInverseBoundaryMap (int eid)
virtual Array< OneD, unsigned int > v_GetFaceInverseBoundaryMap (int fid, StdRegions::Orientation faceOrient=StdRegions::eNoOrientation)
virtual DNekMatSharedPtr v_BuildTransformationMatrix (const DNekScalMatSharedPtr &r_bnd, const StdRegions::MatrixType matrixType)
virtual DNekMatSharedPtr v_BuildInverseTransformationMatrix (const DNekScalMatSharedPtr &m_transformationmatrix)
 Build inverse and inverse transposed transformation matrix: $\mathbf{R^{-1}}$ and $\mathbf{R^{-T}}$.
virtual DNekMatSharedPtr v_BuildVertexMatrix (const DNekScalMatSharedPtr &r_bnd)
- Protected Member Functions inherited from Nektar::LocalRegions::Expansion
void ComputeLaplacianMetric ()
void ComputeQuadratureMetric ()
virtual void v_MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_ComputeLaplacianMetric ()

Private Member Functions

virtual void v_LaplacianMatrixOp_MatFree_Kernel (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
 Calculate the Laplacian multiplication in a matrix-free manner.

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::StdExpansion3D
std::map< int, NormalVectorm_faceNormals
std::map< int, bool > m_negatedNormals
- Protected Attributes inherited from Nektar::LocalRegions::Expansion
SpatialDomains::GeometrySharedPtr m_geom
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
MetricMap m_metrics

Detailed Description

Definition at line 50 of file PrismExp.h.

Constructor & Destructor Documentation

Nektar::LocalRegions::PrismExp::PrismExp ( const LibUtilities::BasisKey Ba,
const LibUtilities::BasisKey Bb,
const LibUtilities::BasisKey Bc,
const SpatialDomains::PrismGeomSharedPtr geom 
)

Constructor using BasisKey class for quadrature points and order definition.

Definition at line 47 of file PrismExp.cpp.

:
Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
3, Ba, Bb, Bc),
Ba.GetNumModes(), Bb.GetNumModes(), Bc.GetNumModes()),
Ba, Bb, Bc),
StdPrismExp (Ba, Bb, Bc),
Expansion (geom),
Expansion3D (geom),
boost::bind(&PrismExp::CreateMatrix, this, _1),
std::string("PrismExpMatrix")),
boost::bind(&PrismExp::CreateStaticCondMatrix, this, _1),
std::string("PrismExpStaticCondMatrix"))
{
}
Nektar::LocalRegions::PrismExp::PrismExp ( const PrismExp T)

Definition at line 69 of file PrismExp.cpp.

:
m_matrixManager(T.m_matrixManager),
m_staticCondMatrixManager(T.m_staticCondMatrixManager)
{
}
Nektar::LocalRegions::PrismExp::~PrismExp ( )

Definition at line 80 of file PrismExp.cpp.

{
}

Member Function Documentation

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

Definition at line 1308 of file PrismExp.cpp.

References ASSERTL2, Nektar::LocalRegions::Expansion::BuildTransformationMatrix(), Nektar::LocalRegions::Expansion::BuildVertexMatrix(), Nektar::StdRegions::StdMatrixKey::ConstFactorExists(), Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::eHelmholtz, Nektar::StdRegions::eHybridDGHelmBndLam, Nektar::StdRegions::eHybridDGHelmholtz, Nektar::StdRegions::eHybridDGLamToQ0, Nektar::StdRegions::eHybridDGLamToQ1, Nektar::StdRegions::eHybridDGLamToU, Nektar::StdRegions::eInvHybridDGHelmholtz, Nektar::StdRegions::eInvLaplacianWithUnityMean, Nektar::StdRegions::eInvMass, Nektar::StdRegions::eIProductWRTBase, Nektar::StdRegions::eLaplacian, Nektar::StdRegions::eLaplacian00, Nektar::StdRegions::eLaplacian01, Nektar::StdRegions::eLaplacian02, Nektar::StdRegions::eLaplacian11, Nektar::StdRegions::eLaplacian12, Nektar::StdRegions::eLaplacian22, Nektar::StdRegions::eMass, Nektar::SpatialDomains::eNoGeomType, Nektar::StdRegions::ePreconLinearSpace, Nektar::StdRegions::ePreconLinearSpaceMass, Nektar::StdRegions::ePreconR, Nektar::StdRegions::ePreconRMass, Nektar::StdRegions::ePreconRT, Nektar::StdRegions::ePreconRTMass, 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(), m_matrixManager, Nektar::LocalRegions::Expansion::m_metricinfo, and Nektar::Transpose().

{
"Geometric information is not set up");
switch(mkey.GetMatrixType())
{
{
{
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;
StdRegions::StdMatrixKey masskey(StdRegions::eMass,DetShapeType(),*this);
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;
}
{
{
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);
MatrixKey deriv2key(StdRegions::eWeakDeriv2,
mkey.GetShapeType(), *this);
DNekMat &deriv0 = *GetStdMatrix(deriv0key);
DNekMat &deriv1 = *GetStdMatrix(deriv1key);
DNekMat &deriv2 = *GetStdMatrix(deriv2key);
int rows = deriv0.GetRows();
int cols = deriv1.GetColumns();
DNekMatSharedPtr WeakDeriv = MemoryManager<DNekMat>
::AllocateSharedPtr(rows,cols);
(*WeakDeriv) = df[3*dir ][0]*deriv0
+ df[3*dir+1][0]*deriv1
+ df[3*dir+2][0]*deriv2;
returnval = MemoryManager<DNekScalMat>
::AllocateSharedPtr(jac,WeakDeriv);
}
break;
}
{
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 lap02key(StdRegions::eLaplacian02,
mkey.GetShapeType(), *this);
MatrixKey lap11key(StdRegions::eLaplacian11,
mkey.GetShapeType(), *this);
MatrixKey lap12key(StdRegions::eLaplacian12,
mkey.GetShapeType(), *this);
MatrixKey lap22key(StdRegions::eLaplacian22,
mkey.GetShapeType(), *this);
DNekMat &lap00 = *GetStdMatrix(lap00key);
DNekMat &lap01 = *GetStdMatrix(lap01key);
DNekMat &lap02 = *GetStdMatrix(lap02key);
DNekMat &lap11 = *GetStdMatrix(lap11key);
DNekMat &lap12 = *GetStdMatrix(lap12key);
DNekMat &lap22 = *GetStdMatrix(lap22key);
NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
Array<TwoD, const NekDouble> gmat
= m_metricinfo->GetGmat(ptsKeys);
int rows = lap00.GetRows();
int cols = lap00.GetColumns();
DNekMatSharedPtr lap = MemoryManager<DNekMat>
::AllocateSharedPtr(rows,cols);
(*lap) = gmat[0][0]*lap00
+ gmat[4][0]*lap11
+ gmat[8][0]*lap22
+ gmat[3][0]*(lap01 + Transpose(lap01))
+ gmat[6][0]*(lap02 + Transpose(lap02))
+ gmat[7][0]*(lap12 + Transpose(lap12));
returnval = MemoryManager<DNekScalMat>
::AllocateSharedPtr(jac,lap);
}
break;
}
{
NekDouble factor = mkey.GetConstFactor(StdRegions::eFactorLambda);
MatrixKey masskey(StdRegions::eMass,
mkey.GetShapeType(), *this);
DNekScalMat &MassMat = *(this->m_matrixManager[masskey]);
MatrixKey lapkey(StdRegions::eLaplacian,
mkey.GetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
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 + factor*MassMat;
returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,helm);
break;
}
{
NekDouble one = 1.0;
returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,mat);
break;
}
{
NekDouble one = 1.0;
MatrixKey hkey(StdRegions::eHybridDGHelmholtz, DetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
// StdRegions::StdMatrixKey hkey(StdRegions::eHybridDGHelmholtz,
// DetExpansionType(),*this,
// mkey.GetConstant(0),
// mkey.GetConstant(1));
mat->Invert();
returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,mat);
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;
MatrixKey helmkey(StdRegions::eHelmholtz, mkey.GetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
DNekScalMatSharedPtr A =helmStatCond->GetBlock(0,0);
returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,R);
}
break;
{
NekDouble one = 1.0;
MatrixKey masskey(StdRegions::eMass, mkey.GetShapeType(), *this);
DNekScalMatSharedPtr A =massStatCond->GetBlock(0,0);
returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,R);
}
break;
{
NekDouble one = 1.0;
MatrixKey helmkey(StdRegions::eHelmholtz, mkey.GetShapeType(), *this,mkey.GetConstFactors(), mkey.GetVarCoeffs());
DNekScalMatSharedPtr A =helmStatCond->GetBlock(0,0);
DNekMatSharedPtr R=BuildTransformationMatrix(A,mkey.GetMatrixType());
returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,R);
}
break;
{
NekDouble one = 1.0;
MatrixKey helmkey(StdRegions::eHelmholtz, mkey.GetShapeType(), *this,mkey.GetConstFactors(), mkey.GetVarCoeffs());
DNekScalMatSharedPtr A =helmStatCond->GetBlock(0,0);
DNekMatSharedPtr R=BuildTransformationMatrix(A,mkey.GetMatrixType());
returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,R);
}
break;
{
NekDouble one = 1.0;
MatrixKey masskey(StdRegions::eMass, mkey.GetShapeType(), *this);
DNekScalMatSharedPtr A =massStatCond->GetBlock(0,0);
DNekMatSharedPtr R=BuildTransformationMatrix(A,mkey.GetMatrixType());
returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,R);
}
break;
{
NekDouble one = 1.0;
MatrixKey masskey(StdRegions::eMass, mkey.GetShapeType(), *this);
DNekScalMatSharedPtr A =massStatCond->GetBlock(0,0);
DNekMatSharedPtr R=BuildTransformationMatrix(A,mkey.GetMatrixType());
returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,R);
}
break;
default:
{
NekDouble one = 1.0;
returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,mat);
}
}
return returnval;
}
DNekScalBlkMatSharedPtr Nektar::LocalRegions::PrismExp::CreateStaticCondMatrix ( const MatrixKey mkey)
protected

Definition at line 1625 of file PrismExp.cpp.

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

{
ASSERTL2(m_metricinfo->GetGtype() != SpatialDomains::eNoGeomType,"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);
NekDouble factor = 1.0;
switch(mkey.GetMatrixType())
{
case StdRegions::eHelmholtz: // special case since Helmholtz not defined in StdRegions
// use Deformed case for both regular and deformed geometries
factor = 1.0;
goto UseLocRegionsMatrix;
break;
default:
{
factor = 1.0;
goto UseLocRegionsMatrix;
}
else
{
factor = mat->Scale();
goto UseStdRegionsMatrix;
}
break;
UseStdRegionsMatrix:
{
NekDouble invfactor = 1.0/factor;
NekDouble one = 1.0;
//TODO: check below
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));
}
break;
}
return returnval;
}
void Nektar::LocalRegions::PrismExp::v_ComputeFaceNormal ( const int  face)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 869 of file PrismExp.cpp.

References ASSERTL0, Nektar::SpatialDomains::eMovingRegular, Nektar::SpatialDomains::eRegular, Vmath::Fill(), Nektar::StdRegions::StdExpansion::GetCoordim(), Nektar::LocalRegions::Expansion::GetGeom(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::LibUtilities::Interp2D(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion3D::m_faceNormals, Vmath::Sdiv(), Vmath::Smul(), Vmath::Vmul(), Vmath::Vsqrt(), Vmath::Vvtvp(), and Vmath::Zero().

{
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);
// Number of quadrature points in face expansion.
int nq = m_base[0]->GetNumPoints()*m_base[0]->GetNumPoints();
int vCoordDim = GetCoordim();
int i;
m_faceNormals[face] = Array<OneD, Array<OneD, NekDouble> >(vCoordDim);
Array<OneD, Array<OneD, NekDouble> > &normal = m_faceNormals[face];
for (i = 0; i < vCoordDim; ++i)
{
normal[i] = Array<OneD, NekDouble>(nq);
}
// Regular geometry case
if (type == SpatialDomains::eRegular ||
{
NekDouble fac;
// Set up normals
switch(face)
{
case 0:
{
for(i = 0; i < vCoordDim; ++i)
{
Vmath::Fill(nq,-df[3*i+2][0],normal[i],1);
}
break;
}
case 1:
{
for(i = 0; i < vCoordDim; ++i)
{
Vmath::Fill(nq,-df[3*i+1][0],normal[i],1);
}
break;
}
case 2:
{
for(i = 0; i < vCoordDim; ++i)
{
Vmath::Fill(nq,df[3*i][0]+df[3*i+2][0],normal[i],1);
}
break;
}
case 3:
{
for(i = 0; i < vCoordDim; ++i)
{
Vmath::Fill(nq,df[3*i+1][0],normal[i],1);
}
break;
}
case 4:
{
for(i = 0; i < vCoordDim; ++i)
{
Vmath::Fill(nq,-df[3*i][0],normal[i],1);
}
break;
}
default:
ASSERTL0(false,"face is out of range (face < 4)");
}
// Normalise resulting vector.
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(nq,fac,normal[i],1,normal[i],1);
}
}
else
{
// Set up deformed normals.
int j, k;
int nq0 = ptsKeys[0].GetNumPoints();
int nq1 = ptsKeys[1].GetNumPoints();
int nq2 = ptsKeys[2].GetNumPoints();
int nq01 = nq0*nq1;
int nqtot;
// Determine number of quadrature points on the face.
if (face == 0)
{
nqtot = nq0*nq1;
}
else if (face == 1 || face == 3)
{
nqtot = nq0*nq2;
}
else
{
nqtot = nq1*nq2;
}
LibUtilities::PointsKey points0;
LibUtilities::PointsKey points1;
Array<OneD, NekDouble> work (nq, 0.0);
Array<OneD, NekDouble> normals(vCoordDim*nqtot,0.0);
// Extract Jacobian along face and recover local derivatives
// (dx/dr) for polynomial interpolation by multiplying m_gmat by
// jacobian
switch(face)
{
case 0:
{
for(j = 0; j < nq01; ++j)
{
normals[j] = -df[2][j]*jac[j];
normals[nqtot+j] = -df[5][j]*jac[j];
normals[2*nqtot+j] = -df[8][j]*jac[j];
}
points0 = ptsKeys[0];
points1 = ptsKeys[1];
break;
}
case 1:
{
for (j = 0; j < nq0; ++j)
{
for(k = 0; k < nq2; ++k)
{
int tmp = j+nq01*k;
normals[j+k*nq0] =
-df[1][tmp]*jac[tmp];
normals[nqtot+j+k*nq0] =
-df[4][tmp]*jac[tmp];
normals[2*nqtot+j+k*nq0] =
-df[7][tmp]*jac[tmp];
}
}
points0 = ptsKeys[0];
points1 = ptsKeys[2];
break;
}
case 2:
{
for (j = 0; j < nq1; ++j)
{
for(k = 0; k < nq2; ++k)
{
int tmp = nq0-1+nq0*j+nq01*k;
normals[j+k*nq1] =
(df[0][tmp]+df[2][tmp])*jac[tmp];
normals[nqtot+j+k*nq1] =
(df[3][tmp]+df[5][tmp])*jac[tmp];
normals[2*nqtot+j+k*nq1] =
(df[6][tmp]+df[8][tmp])*jac[tmp];
}
}
points0 = ptsKeys[1];
points1 = ptsKeys[2];
break;
}
case 3:
{
for (j = 0; j < nq0; ++j)
{
for(k = 0; k < nq2; ++k)
{
int tmp = nq0*(nq1-1) + j + nq01*k;
normals[j+k*nq0] =
df[1][tmp]*jac[tmp];
normals[nqtot+j+k*nq0] =
df[4][tmp]*jac[tmp];
normals[2*nqtot+j+k*nq0] =
df[7][tmp]*jac[tmp];
}
}
points0 = ptsKeys[0];
points1 = ptsKeys[2];
break;
}
case 4:
{
for (j = 0; j < nq1; ++j)
{
for(k = 0; k < nq2; ++k)
{
int tmp = j*nq0+nq01*k;
normals[j+k*nq1] =
-df[0][tmp]*jac[tmp];
normals[nqtot+j+k*nq1] =
-df[3][tmp]*jac[tmp];
normals[2*nqtot+j+k*nq1] =
-df[6][tmp]*jac[tmp];
}
}
points0 = ptsKeys[1];
points1 = ptsKeys[2];
break;
}
default:
ASSERTL0(false,"face is out of range (face < 4)");
}
// Interpolate Jacobian and invert
LibUtilities::Interp2D(points0, points1, jac,
m_base[0]->GetPointsKey(),
m_base[0]->GetPointsKey(),
work);
Vmath::Sdiv(nq, 1.0, &work[0], 1, &work[0], 1);
// Interpolate normal and multiply by inverse Jacobian.
for(i = 0; i < vCoordDim; ++i)
{
LibUtilities::Interp2D(points0, points1,
&normals[i*nqtot],
m_base[0]->GetPointsKey(),
m_base[0]->GetPointsKey(),
&normal[i][0]);
Vmath::Vmul(nq,work,1,normal[i],1,normal[i],1);
}
// Normalise to obtain unit normals.
Vmath::Zero(nq,work,1);
for(i = 0; i < GetCoordim(); ++i)
{
Vmath::Vvtvp(nq,normal[i],1,normal[i],1,work,1,work,1);
}
Vmath::Vsqrt(nq,work,1,work,1);
Vmath::Sdiv (nq,1.0,work,1,work,1);
for(i = 0; i < GetCoordim(); ++i)
{
Vmath::Vmul(nq,normal[i],1,work,1,normal[i],1);
}
}
}
DNekMatSharedPtr Nektar::LocalRegions::PrismExp::v_CreateStdMatrix ( const StdRegions::StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdPrismExp.

Definition at line 1282 of file PrismExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1303 of file PrismExp.cpp.

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

void Nektar::LocalRegions::PrismExp::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 501 of file PrismExp.cpp.

References ASSERTL0, ASSERTL1, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Vmath::Vcopy(), and Vmath::Zero().

{
int data_order0 = nummodes[mode_offset];
int fillorder0 = 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);
int data_order2 = nummodes[mode_offset+2];
int order2 = m_base[2]->GetNumModes();
int fillorder2 = min(order2,data_order2);
switch(m_base[0]->GetBasisType())
{
{
int i,j;
int cnt = 0;
int cnt1 = 0;
"Extraction routine not set up for this basis");
"Extraction routine not set up for this basis");
for(j = 0; j < fillorder0; ++j)
{
for(i = 0; i < fillorder1; ++i)
{
Vmath::Vcopy(fillorder2-j, &data[cnt], 1,
&coeffs[cnt1], 1);
cnt += data_order2-j;
cnt1 += order2-j;
}
// count out data for j iteration
for(i = fillorder1; i < data_order1; ++i)
{
cnt += data_order2-j;
}
for(i = fillorder1; i < order1; ++i)
{
cnt1 += order2-j;
}
}
}
break;
default:
ASSERTL0(false, "basis is either not set up or not "
"hierarchicial");
}
}
void Nektar::LocalRegions::PrismExp::v_FwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Forward transform from physical quadrature space stored in inarray and evaluate the expansion coefficients and store in (this)->m_coeffs.

Inputs:

Outputs:

Reimplemented from Nektar::StdRegions::StdPrismExp.

Definition at line 215 of file PrismExp.cpp.

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

{
if(m_base[0]->Collocation() &&
m_base[1]->Collocation() &&
m_base[2]->Collocation())
{
Vmath::Vcopy(GetNcoeffs(),&inarray[0],1,&outarray[0],1);
}
else
{
v_IProductWRTBase(inarray, outarray);
// get Mass matrix inverse
MatrixKey masskey(StdRegions::eInvMass,
DetShapeType(),*this);
// copy inarray in case inarray == outarray
DNekVec in (m_ncoeffs,outarray);
DNekVec out(m_ncoeffs,outarray,eWrapper);
out = (*matsys)*in;
}
}
void Nektar::LocalRegions::PrismExp::v_GeneralMatrixOp_MatOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdRegions::StdMatrixKey mkey 
)
protectedvirtual

Definition at line 1161 of file PrismExp.cpp.

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

{
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::PrismExp::v_GenMatrix ( const StdRegions::StdMatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdPrismExp.

Definition at line 1259 of file PrismExp.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().

void Nektar::LocalRegions::PrismExp::v_GetCoord ( const Array< OneD, const NekDouble > &  Lcoords,
Array< OneD, NekDouble > &  coords 
)
protectedvirtual

Get the coordinates #coords at the local coordinates #Lcoords.

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 440 of file PrismExp.cpp.

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

{
int i;
ASSERTL1(Lcoords[0] <= -1.0 && Lcoords[0] >= 1.0 &&
Lcoords[1] <= -1.0 && Lcoords[1] >= 1.0 &&
Lcoords[2] <= -1.0 && Lcoords[2] >= 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::PrismExp::v_GetCoordim ( void  )
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 496 of file PrismExp.cpp.

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

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

Reimplemented from Nektar::StdRegions::StdPrismExp.

Definition at line 458 of file PrismExp.cpp.

{
Expansion::v_GetCoords(coords_0, coords_1, coords_2);
}
StdRegions::Orientation Nektar::LocalRegions::PrismExp::v_GetFaceOrient ( int  face)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 562 of file PrismExp.cpp.

References Nektar::LocalRegions::Expansion3D::GetGeom3D().

{
return GetGeom3D()->GetFaceOrient(face);
}
void Nektar::LocalRegions::PrismExp::v_GetFacePhysVals ( const int  face,
const StdRegions::StdExpansionSharedPtr FaceExp,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
StdRegions::Orientation  orient 
)
protectedvirtual

Definition at line 578 of file PrismExp.cpp.

References ASSERTL0, Nektar::StdRegions::eDir1BwdDir1_Dir2BwdDir2, Nektar::StdRegions::eDir1BwdDir1_Dir2FwdDir2, Nektar::StdRegions::eDir1BwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1BwdDir2_Dir2FwdDir1, Nektar::StdRegions::eDir1FwdDir1_Dir2BwdDir2, Nektar::StdRegions::eDir1FwdDir1_Dir2FwdDir2, Nektar::StdRegions::eDir1FwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1FwdDir2_Dir2FwdDir1, Nektar::StdRegions::eNoOrientation, Nektar::StdRegions::StdExpansion::GetFaceOrient(), Nektar::LibUtilities::Interp2D(), Nektar::StdRegions::StdExpansion::m_base, and Vmath::Vcopy().

Referenced by v_GetTracePhysVals().

{
int nquad0 = m_base[0]->GetNumPoints();
int nquad1 = m_base[1]->GetNumPoints();
int nquad2 = m_base[2]->GetNumPoints();
Array<OneD,NekDouble> o_tmp(nquad0*nquad1*nquad2);
{
orient = GetFaceOrient(face);
}
switch(face)
{
case 0:
{
//Directions A and B positive
Vmath::Vcopy(nquad0*nquad1,&(inarray[0]),1,&(outarray[0]),1);
}
{
//Direction A negative and B positive
for (int j=0; j<nquad1; j++)
{
Vmath::Vcopy(nquad0,&(inarray[0])+(nquad0-1)+j*nquad0,-1,&(outarray[0])+(j*nquad0),1);
}
}
{
//Direction A positive and B negative
for (int j=0; j<nquad1; j++)
{
Vmath::Vcopy(nquad0,&(inarray[0])+nquad0*(nquad1-1-j),1,&(outarray[0])+(j*nquad0),1);
}
}
{
//Direction A negative and B negative
for(int j=0; j<nquad1; j++)
{
Vmath::Vcopy(nquad0,&(inarray[0])+(nquad0*nquad1-1-j*nquad0),-1,&(outarray[0])+(j*nquad0),1);
}
}
{
//Transposed, Direction A and B positive
for (int i=0; i<nquad0; i++)
{
Vmath::Vcopy(nquad1,&(inarray[0])+i,nquad0,&(outarray[0])+(i*nquad1),1);
}
}
{
//Transposed, Direction A positive and B negative
for (int i=0; i<nquad0; i++)
{
Vmath::Vcopy(nquad1,&(inarray[0])+(nquad0-1-i),nquad0,&(outarray[0])+(i*nquad1),1);
}
}
{
//Transposed, Direction A negative and B positive
for (int i=0; i<nquad0; i++)
{
Vmath::Vcopy(nquad1,&(inarray[0])+i+nquad0*(nquad1-1),-nquad0,&(outarray[0])+(i*nquad1),1);
}
}
{
//Transposed, Direction A and B negative
for (int i=0; i<nquad0; i++)
{
Vmath::Vcopy(nquad1,&(inarray[0])+(nquad0*nquad1-1-i),-nquad0,&(outarray[0])+(i*nquad1),1);
}
}
o_tmp=outarray;
//interpolate
LibUtilities::Interp2D(m_base[0]->GetPointsKey(), m_base[1]->GetPointsKey(), o_tmp,
FaceExp->GetBasis(0)->GetPointsKey(),FaceExp->GetBasis(1)->GetPointsKey(),outarray);
break;
case 1:
{
//Direction A and B positive
for (int k=0; k<nquad2; k++)
{
Vmath::Vcopy(nquad0,&(inarray[0])+(nquad0*nquad1*k),1,&(outarray[0])+(k*nquad0),1);
}
o_tmp=outarray;
}
else
{
//Direction A negative and B positive
for (int k=0; k<nquad2; k++)
{
Vmath::Vcopy(nquad0,&(inarray[0])+(nquad0-1)+(nquad0*nquad1*k),-1,&(outarray[0])+(k*nquad0),1);
}
o_tmp=outarray;
}
//interpolate
LibUtilities::Interp2D(m_base[0]->GetPointsKey(), m_base[2]->GetPointsKey(), o_tmp,
FaceExp->GetBasis(0)->GetPointsKey(),FaceExp->GetBasis(1)->GetPointsKey(),outarray);
break;
case 2:
{
//Directions A and B positive
Vmath::Vcopy(nquad0*nquad2,&(inarray[0])+(nquad0-1),
nquad0,&(outarray[0]),1);
}
{
//Direction A negative and B positive
for (int k=0; k<nquad2; k++)
{
Vmath::Vcopy(nquad0,&(inarray[0])+(nquad0*nquad1-1)+(k*nquad0*nquad1),
-nquad0,&(outarray[0])+(k*nquad0),1);
}
}
{
//Direction A positive and B negative
for (int k=0; k<nquad2; k++)
{
Vmath::Vcopy(nquad0,&(inarray[0])+(nquad0-1)+(nquad0*nquad1*(nquad2-1-k)),
nquad0,&(outarray[0])+(k*nquad0),1);
}
}
{
//Direction A negative and B negative
for (int k=0; k<nquad2; k++)
{
Vmath::Vcopy(nquad0,&(inarray[0])+(nquad0*nquad1-1)+(nquad0*nquad1*(nquad2-1-k)),
-nquad0,&(outarray[0])+(k*nquad0),1);
}
}
{
//Transposed, Direction A and B positive
for (int j=0; j<nquad1; j++)
{
Vmath::Vcopy(nquad2,&(inarray[0])+(nquad0-1)+(j*nquad0),
nquad0*nquad1,&(outarray[0])+(j*nquad2),1);
}
}
{
//Transposed, Direction A positive and B negative
for (int j=0; j<nquad0; j++)
{
Vmath::Vcopy(nquad2,&(inarray[0])+(nquad0*nquad1-1-j*nquad0),
nquad0*nquad1,&(outarray[0])+(j*nquad2),1);
}
}
{
//Transposed, Direction A negative and B positive
for (int j=0; j<nquad0; j++)
{
Vmath::Vcopy(nquad2,&(inarray[0])+nquad0*nquad1*(nquad2-1)+nquad0+j*nquad0,
-nquad0*nquad1,&(outarray[0])+(j*nquad2),1);
}
}
{
//Transposed, Direction A and B negative
for (int j=0; j<nquad0; j++)
{
Vmath::Vcopy(nquad2,&(inarray[0])+(nquad0*nquad1*nquad2-1-j*nquad0),
-nquad0*nquad1,&(outarray[0])+(j*nquad2),1);
}
}
o_tmp=outarray;
//interpolate
LibUtilities::Interp2D(m_base[1]->GetPointsKey(), m_base[2]->GetPointsKey(), o_tmp,
FaceExp->GetBasis(0)->GetPointsKey(),FaceExp->GetBasis(1)->GetPointsKey(),outarray);
break;
case 3:
{
//Directions A and B positive
for (int k=0; k<nquad2; k++)
{
Vmath::Vcopy(nquad0,&(inarray[0])+(nquad0*(nquad1-1))+(k*nquad0*nquad1),
1,&(outarray[0])+(k*nquad0),1);
}
}
else
{
//Direction A negative and B positive
for (int k=0; k<nquad2; k++)
{
Vmath::Vcopy(nquad0,&(inarray[0])+(nquad0*nquad1-1)+(k*nquad0*nquad1),
-1,&(outarray[0])+(k*nquad0),1);
}
}
o_tmp=outarray;
//interpolate
LibUtilities::Interp2D(m_base[0]->GetPointsKey(), m_base[2]->GetPointsKey(), o_tmp,
FaceExp->GetBasis(0)->GetPointsKey(),FaceExp->GetBasis(1)->GetPointsKey(),outarray);
break;
case 4:
{
//Directions A and B positive
Vmath::Vcopy(nquad1*nquad2,&(inarray[0]),nquad0,&(outarray[0]),1);
}
{
//Direction A negative and B positive
for (int k=0; k<nquad2; k++)
{
Vmath::Vcopy(nquad1,&(inarray[0])+nquad0*(nquad1-1)+(k*nquad0*nquad1),
-nquad0,&(outarray[0])+(k*nquad1),1);
}
}
{
//Direction A positive and B negative
for (int k=0; k<nquad2; k++)
{
Vmath::Vcopy(nquad1,&(inarray[0])+(nquad0*nquad1*(nquad2-1-k)),
nquad0,&(outarray[0])+(k*nquad1),1);
}
}
{
//Direction A negative and B negative
for (int k=0; k<nquad2; k++)
{
Vmath::Vcopy(nquad1,&(inarray[0])+nquad0*(nquad1-1)+(nquad0*nquad1*(nquad2-1-k)),
-nquad0,&(outarray[0])+(k*nquad1),1);
}
}
{
//Transposed, Direction A and B positive
for (int j=0; j<nquad1; j++)
{
Vmath::Vcopy(nquad2,&(inarray[0])+j*nquad0,nquad0*nquad1,
&(outarray[0])+(j*nquad2),1);
}
}
{
//Transposed, Direction A positive and B negative
for (int j=0; j<nquad1; j++)
{
Vmath::Vcopy(nquad2,&(inarray[0])+(nquad0*(nquad1-1)-j*nquad0),
nquad0*nquad1,&(outarray[0])+(j*nquad2),1);
}
}
{
//Transposed, Direction A negative and B positive
for (int j=0; j<nquad1; j++)
{
Vmath::Vcopy(nquad2,&(inarray[0])+nquad0*nquad1*(nquad2-1)+j*nquad0,
-nquad0*nquad1,&(outarray[0])+(j*nquad2),1);
}
}
{
//Transposed, Direction A and B negative
for (int j=0; j<nquad1; j++)
{
Vmath::Vcopy(nquad2,&(inarray[0])+(nquad0*(nquad1*nquad2-1)-j*nquad0),
-nquad0*nquad1,&(outarray[0])+(j*nquad2),1);
}
}
o_tmp=outarray;
//interpolate
LibUtilities::Interp2D(m_base[1]->GetPointsKey(), m_base[2]->GetPointsKey(), o_tmp,
FaceExp->GetBasis(0)->GetPointsKey(),FaceExp->GetBasis(1)->GetPointsKey(),outarray);
break;
default:
ASSERTL0(false,"face value (> 4) is out of range");
break;
}
}
DNekScalMatSharedPtr Nektar::LocalRegions::PrismExp::v_GetLocMatrix ( const MatrixKey mkey)
protectedvirtual

Reimplemented from Nektar::LocalRegions::Expansion.

Definition at line 1293 of file PrismExp.cpp.

References m_matrixManager.

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1298 of file PrismExp.cpp.

References m_staticCondMatrixManager.

{
}
void Nektar::LocalRegions::PrismExp::v_GetTracePhysVals ( const int  face,
const StdRegions::StdExpansionSharedPtr FaceExp,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
StdRegions::Orientation  orient 
)
protectedvirtual

Returns the physical values at the quadrature points of a face.

Definition at line 568 of file PrismExp.cpp.

References v_GetFacePhysVals().

{
v_GetFacePhysVals(face,FaceExp,inarray,outarray,orient);
}
void Nektar::LocalRegions::PrismExp::v_HelmholtzMatrixOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdRegions::StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1153 of file PrismExp.cpp.

References Nektar::StdRegions::StdExpansion3D::v_HelmholtzMatrixOp_MatFree().

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

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

Inputs:

Outputs:

Reimplemented from Nektar::StdRegions::StdExpansion3D.

Definition at line 109 of file PrismExp.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();
int nquad2 = m_base[2]->GetNumPoints();
Array<OneD, const NekDouble> jac = m_metricinfo->GetJac(GetPointsKeys());
Array<OneD, NekDouble> tmp(nquad0*nquad1*nquad2);
// Multiply inarray with Jacobian
{
Vmath::Vmul(nquad0*nquad1*nquad2,&jac[0],1,(NekDouble*)&inarray[0],1,&tmp[0],1);
}
else
{
Vmath::Smul(nquad0*nquad1*nquad2,(NekDouble)jac[0],(NekDouble*)&inarray[0],1,&tmp[0],1);
}
// Call StdPrismExp version.
}
void Nektar::LocalRegions::PrismExp::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*base2 and put into outarray:

$ \begin{array}{rcl} I_{pqr} = (\phi_{pqr}, u)_{\delta} & = & \sum_{i=0}^{nq_0} \sum_{j=0}^{nq_1} \sum_{k=0}^{nq_2} \psi_{p}^{a} (\bar \eta_{1i}) \psi_{q}^{a} (\xi_{2j}) \psi_{pr}^{b} (\xi_{3k}) w_i w_j w_k u(\bar \eta_{1,i} \xi_{2,j} \xi_{3,k}) J_{i,j,k}\\ & = & \sum_{i=0}^{nq_0} \psi_p^a(\bar \eta_{1,i}) \sum_{j=0}^{nq_1} \psi_{q}^a(\xi_{2,j}) \sum_{k=0}^{nq_2} \psi_{pr}^b u(\bar \eta_{1i},\xi_{2j},\xi_{3k}) J_{i,j,k} \end{array} $

where

$ \phi_{pqr} (\xi_1 , \xi_2 , \xi_3) = \psi_p^a (\bar \eta_1) \psi_{q}^a (\xi_2) \psi_{pr}^b (\xi_3) $

which can be implemented as
$f_{pr} (\xi_{3k}) = \sum_{k=0}^{nq_3} \psi_{pr}^b u(\bar \eta_{1i},\xi_{2j},\xi_{3k}) J_{i,j,k} = {\bf B_3 U} $
$ g_{q} (\xi_{3k}) = \sum_{j=0}^{nq_1} \psi_{q}^a (\xi_{2j}) f_{pr} (\xi_{3k}) = {\bf B_2 F} $
$ (\phi_{pqr}, u)_{\delta} = \sum_{k=0}^{nq_0} \psi_{p}^a (\xi_{3k}) g_{q} (\xi_{3k}) = {\bf B_1 G} $

Reimplemented from Nektar::StdRegions::StdPrismExp.

Definition at line 270 of file PrismExp.cpp.

References v_IProductWRTBase_SumFac().

Referenced by v_FwdTrans().

{
v_IProductWRTBase_SumFac(inarray, outarray);
}
void Nektar::LocalRegions::PrismExp::v_IProductWRTBase_SumFac ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdPrismExp.

Definition at line 277 of file PrismExp.cpp.

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

Referenced by v_IProductWRTBase().

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

Calculates the inner product $ I_{pqr} = (u, \partial_{x_i} \phi_{pqr}) $.

The derivative of the basis functions is performed using the chain rule in order to incorporate the geometric factors. Assuming that the basis functions are a tensor product $\phi_{pqr}(\eta_1,\eta_2,\eta_3) = \phi_1(\eta_1)\phi_2(\eta_2)\phi_3(\eta_3)$, this yields the result

\[ I_{pqr} = \sum_{j=1}^3 \left(u, \frac{\partial u}{\partial \eta_j} \frac{\partial \eta_j}{\partial x_i}\right) \]

In the tetrahedral element, we must also incorporate a second set of geometric factors which incorporate the collapsed co-ordinate system, so that

\[ \frac{\partial\eta_j}{\partial x_i} = \sum_{k=1}^3 \frac{\partial\eta_j}{\partial\xi_k}\frac{\partial\xi_k}{\partial x_i} \]

These derivatives can be found on p152 of Sherwin & Karniadakis.

Parameters
dirDirection in which to take the derivative.
inarrayThe function $ u $.
outarrayValue of the inner product.

Reimplemented from Nektar::StdRegions::StdPrismExp.

Definition at line 329 of file PrismExp.cpp.

References v_IProductWRTDerivBase_SumFac().

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

Reimplemented from Nektar::StdRegions::StdPrismExp.

Definition at line 337 of file PrismExp.cpp.

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

Referenced by v_IProductWRTDerivBase().

{
const int nquad0 = m_base[0]->GetNumPoints();
const int nquad1 = m_base[1]->GetNumPoints();
const int nquad2 = m_base[2]->GetNumPoints();
const int order0 = m_base[0]->GetNumModes ();
const int order1 = m_base[1]->GetNumModes ();
const int nqtot = nquad0*nquad1*nquad2;
int i;
const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
const Array<OneD, const NekDouble> &z2 = m_base[2]->GetZ();
Array<OneD, NekDouble> gfac0(nquad0 );
Array<OneD, NekDouble> gfac2(nquad2 );
Array<OneD, NekDouble> tmp1 (nqtot );
Array<OneD, NekDouble> tmp2 (nqtot );
Array<OneD, NekDouble> tmp3 (nqtot );
Array<OneD, NekDouble> tmp4 (nqtot );
Array<OneD, NekDouble> tmp5 (nqtot );
Array<OneD, NekDouble> tmp6 (m_ncoeffs);
Array<OneD, NekDouble> wsp (order0*nquad2*(nquad1+order1));
const Array<TwoD, const NekDouble>& df =
m_metricinfo->GetDerivFactors(GetPointsKeys());
{
Vmath::Vmul(nqtot,&df[3*dir][0], 1,tmp1.get(),1,tmp2.get(),1);
Vmath::Vmul(nqtot,&df[3*dir+1][0],1,tmp1.get(),1,tmp3.get(),1);
Vmath::Vmul(nqtot,&df[3*dir+2][0],1,tmp1.get(),1,tmp4.get(),1);
}
else
{
Vmath::Smul(nqtot, df[3*dir][0], tmp1.get(),1,tmp2.get(), 1);
Vmath::Smul(nqtot, df[3*dir+1][0],tmp1.get(),1,tmp3.get(), 1);
Vmath::Smul(nqtot, df[3*dir+2][0],tmp1.get(),1,tmp4.get(), 1);
}
// set up geometric factor: (1+z0)/2
for (i = 0; i < nquad0; ++i)
{
gfac0[i] = 0.5*(1+z0[i]);
}
// Set up geometric factor: 2/(1-z2)
for (i = 0; i < nquad2; ++i)
{
gfac2[i] = 2.0/(1-z2[i]);
}
const int nq01 = nquad0*nquad1;
for (i = 0; i < nquad2; ++i)
{
Vmath::Smul(nq01,gfac2[i],&tmp2[0]+i*nq01,1,&tmp2[0]+i*nq01,1);
Vmath::Smul(nq01,gfac2[i],&tmp4[0]+i*nq01,1,&tmp5[0]+i*nq01,1);
}
for(i = 0; i < nquad1*nquad2; ++i)
{
Vmath::Vmul(nquad0,&gfac0[0],1,&tmp5[0]+i*nquad0,1,
&tmp5[0]+i*nquad0,1);
}
Vmath::Vadd(nqtot, &tmp2[0], 1, &tmp5[0], 1, &tmp2[0], 1);
m_base[1]->GetBdata (),
m_base[2]->GetBdata (),
tmp2,outarray,wsp,
true,true,true);
m_base[1]->GetDbdata(),
m_base[2]->GetBdata (),
tmp3,tmp6,wsp,
true,true,true);
Vmath::Vadd(m_ncoeffs, tmp6, 1, outarray, 1, outarray, 1);
m_base[1]->GetBdata (),
m_base[2]->GetDbdata(),
tmp4,tmp6,wsp,
true,true,true);
Vmath::Vadd(m_ncoeffs, tmp6, 1, outarray, 1, outarray, 1);
}
void Nektar::LocalRegions::PrismExp::v_LaplacianMatrixOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdRegions::StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1135 of file PrismExp.cpp.

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

{
PrismExp::LaplacianMatrixOp_MatFree(inarray,outarray,mkey);
}
void Nektar::LocalRegions::PrismExp::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::StdExpansion.

Definition at line 1143 of file PrismExp.cpp.

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

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

Calculate the Laplacian multiplication in a matrix-free manner.

This function is the kernel of the Laplacian matrix-free operator, and is used in v_HelmholtzMatrixOp_MatFree to determine the effect of the Helmholtz operator in a similar fashion.

The majority of the calculation is precisely the same as in the hexahedral expansion; however the collapsed co-ordinate system must be taken into account when constructing the geometric factors. How this is done is detailed more exactly in the tetrahedral expansion. On entry to this function, the input #inarray must be in its backwards-transformed state (i.e. $\mathbf{u} = \mathbf{B}\hat{\mathbf{u}}$). The output is in coefficient space.

See Also
TetExp::v_HelmholtzMatrixOp_MatFree

Definition at line 1757 of file PrismExp.cpp.

References Nektar::SpatialDomains::eDeformed, Vmath::Fill(), Nektar::StdRegions::StdExpansion::GetPointsKeys(), Nektar::StdRegions::StdExpansion3D::IProductWRTBase_SumFacKernel(), Nektar::StdRegions::StdExpansion::m_base, Nektar::LocalRegions::Expansion::m_metricinfo, Nektar::StdRegions::StdExpansion::m_ncoeffs, Nektar::StdRegions::StdExpansion::MultiplyByQuadratureMetric(), Nektar::StdRegions::StdExpansion3D::PhysTensorDeriv(), Vmath::Svtsvtp(), Vmath::Svtvp(), Vmath::Vadd(), Vmath::Vvtvp(), and Vmath::Vvtvvtp().

{
int nquad0 = m_base[0]->GetNumPoints();
int nquad1 = m_base[1]->GetNumPoints();
int nquad2 = m_base[2]->GetNumPoints();
int nqtot = nquad0*nquad1*nquad2;
int i;
// Set up temporary storage.
Array<OneD,NekDouble> alloc(11*nqtot,0.0);
Array<OneD,NekDouble> wsp1 (alloc ); // TensorDeriv 1
Array<OneD,NekDouble> wsp2 (alloc+ 1*nqtot); // TensorDeriv 2
Array<OneD,NekDouble> wsp3 (alloc+ 2*nqtot); // TensorDeriv 3
Array<OneD,NekDouble> g0 (alloc+ 3*nqtot); // g0
Array<OneD,NekDouble> g1 (alloc+ 4*nqtot); // g1
Array<OneD,NekDouble> g2 (alloc+ 5*nqtot); // g2
Array<OneD,NekDouble> g3 (alloc+ 6*nqtot); // g3
Array<OneD,NekDouble> g4 (alloc+ 7*nqtot); // g4
Array<OneD,NekDouble> g5 (alloc+ 8*nqtot); // g5
Array<OneD,NekDouble> h0 (alloc+ 3*nqtot); // h0 == g0
Array<OneD,NekDouble> h1 (alloc+ 6*nqtot); // h1 == g3
Array<OneD,NekDouble> wsp4 (alloc+ 4*nqtot); // wsp4 == g1
Array<OneD,NekDouble> wsp5 (alloc+ 5*nqtot); // wsp5 == g2
Array<OneD,NekDouble> wsp6 (alloc+ 8*nqtot); // wsp6 == g5
Array<OneD,NekDouble> wsp7 (alloc+ 3*nqtot); // wsp7 == g0
Array<OneD,NekDouble> wsp8 (alloc+ 9*nqtot); // wsp8
Array<OneD,NekDouble> wsp9 (alloc+10*nqtot); // wsp9
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>& base2 = m_base[2]->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>& dbase2 = m_base[2]->GetDbdata();
// Step 1. LAPLACIAN MATRIX OPERATION
// wsp1 = du_dxi1 = D_xi1 * wsp0 = D_xi1 * u
// wsp2 = du_dxi2 = D_xi2 * wsp0 = D_xi2 * u
// wsp3 = du_dxi3 = D_xi3 * wsp0 = D_xi3 * u
StdExpansion3D::PhysTensorDeriv(inarray,wsp1,wsp2,wsp3);
const Array<TwoD, const NekDouble>& df =
m_metricinfo->GetDerivFactors(GetPointsKeys());
const Array<OneD, const NekDouble>& z0 = m_base[0]->GetZ();
const Array<OneD, const NekDouble>& z2 = m_base[2]->GetZ();
// Step 2. Calculate the metric terms of the collapsed
// coordinate transformation (Spencer's book P152)
for (i = 0; i < nquad2; ++i)
{
Vmath::Fill(nquad0*nquad1, 2.0/(1.0-z2[i]), &h0[0]+i*nquad0*nquad1,1);
Vmath::Fill(nquad0*nquad1, 2.0/(1.0-z2[i]), &h1[0]+i*nquad0*nquad1,1);
}
for (i = 0; i < nquad0; i++)
{
Blas::Dscal(nquad1*nquad2, 0.5*(1+z0[i]), &h1[0]+i, nquad0);
}
// Step 3. Construct combined metric terms for physical space to
// collapsed coordinate system. Order of construction optimised
// to minimise temporary storage
{
// wsp4 = d eta_1/d x_1
Vmath::Vvtvvtp(nqtot, &df[0][0], 1, &h0[0], 1, &df[2][0], 1, &h1[0], 1, &wsp4[0], 1);
// wsp5 = d eta_2/d x_1
Vmath::Vvtvvtp(nqtot, &df[3][0], 1, &h0[0], 1, &df[5][0], 1, &h1[0], 1, &wsp5[0], 1);
// wsp6 = d eta_3/d x_1d
Vmath::Vvtvvtp(nqtot, &df[6][0], 1, &h0[0], 1, &df[8][0], 1, &h1[0], 1, &wsp6[0], 1);
// g0 (overwrites h0)
Vmath::Vvtvvtp(nqtot, &wsp4[0], 1, &wsp4[0], 1, &wsp5[0], 1, &wsp5[0], 1, &g0[0], 1);
Vmath::Vvtvp (nqtot, &wsp6[0], 1, &wsp6[0], 1, &g0[0], 1, &g0[0], 1);
// g3 (overwrites h1)
Vmath::Vvtvvtp(nqtot, &df[1][0], 1, &wsp4[0], 1, &df[4][0], 1, &wsp5[0], 1, &g3[0], 1);
Vmath::Vvtvp (nqtot, &df[7][0], 1, &wsp6[0], 1, &g3[0], 1, &g3[0], 1);
// g4
Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &wsp4[0], 1, &df[5][0], 1, &wsp5[0], 1, &g4[0], 1);
Vmath::Vvtvp (nqtot, &df[8][0], 1, &wsp6[0], 1, &g4[0], 1, &g4[0], 1);
// Overwrite wsp4/5/6 with g1/2/5
// g1
Vmath::Vvtvvtp(nqtot, &df[1][0], 1, &df[1][0], 1, &df[4][0], 1, &df[4][0], 1, &g1[0], 1);
Vmath::Vvtvp (nqtot, &df[7][0], 1, &df[7][0], 1, &g1[0], 1, &g1[0], 1);
// g2
Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &df[2][0], 1, &df[5][0], 1, &df[5][0], 1, &g2[0], 1);
Vmath::Vvtvp (nqtot, &df[8][0], 1, &df[8][0], 1, &g2[0], 1, &g2[0], 1);
// g5
Vmath::Vvtvvtp(nqtot, &df[1][0], 1, &df[2][0], 1, &df[4][0], 1, &df[5][0], 1, &g5[0], 1);
Vmath::Vvtvp (nqtot, &df[7][0], 1, &df[8][0], 1, &g5[0], 1, &g5[0], 1);
}
else
{
// wsp4 = d eta_1/d x_1
Vmath::Svtsvtp(nqtot, df[0][0], &h0[0], 1, df[2][0], &h1[0], 1, &wsp4[0], 1);
// wsp5 = d eta_2/d x_1
Vmath::Svtsvtp(nqtot, df[3][0], &h0[0], 1, df[5][0], &h1[0], 1, &wsp5[0], 1);
// wsp6 = d eta_3/d x_1
Vmath::Svtsvtp(nqtot, df[6][0], &h0[0], 1, df[8][0], &h1[0], 1, &wsp6[0], 1);
// g0 (overwrites h0)
Vmath::Vvtvvtp(nqtot, &wsp4[0], 1, &wsp4[0], 1, &wsp5[0], 1, &wsp5[0], 1, &g0[0], 1);
Vmath::Vvtvp (nqtot, &wsp6[0], 1, &wsp6[0], 1, &g0[0], 1, &g0[0], 1);
// g3 (overwrites h1)
Vmath::Svtsvtp(nqtot, df[1][0], &wsp4[0], 1, df[4][0], &wsp5[0], 1, &g3[0], 1);
Vmath::Svtvp (nqtot, df[7][0], &wsp6[0], 1, &g3[0], 1, &g3[0], 1);
// g4
Vmath::Svtsvtp(nqtot, df[2][0], &wsp4[0], 1, df[5][0], &wsp5[0], 1, &g4[0], 1);
Vmath::Svtvp (nqtot, df[8][0], &wsp6[0], 1, &g4[0], 1, &g4[0], 1);
// Overwrite wsp4/5/6 with g1/2/5
// g1
Vmath::Fill(nqtot, df[1][0]*df[1][0] + df[4][0]*df[4][0] + df[7][0]*df[7][0], &g1[0], 1);
// g2
Vmath::Fill(nqtot, df[2][0]*df[2][0] + df[5][0]*df[5][0] + df[8][0]*df[8][0], &g2[0], 1);
// g5
Vmath::Fill(nqtot, df[1][0]*df[2][0] + df[4][0]*df[5][0] + df[7][0]*df[8][0], &g5[0], 1);
}
// Compute component derivatives into wsp7, 8, 9 (wsp7 overwrites
// g0).
Vmath::Vvtvvtp(nqtot,&g0[0],1,&wsp1[0],1,&g3[0],1,&wsp2[0],1,&wsp7[0],1);
Vmath::Vvtvp (nqtot,&g4[0],1,&wsp3[0],1,&wsp7[0],1,&wsp7[0],1);
Vmath::Vvtvvtp(nqtot,&g1[0],1,&wsp2[0],1,&g3[0],1,&wsp1[0],1,&wsp8[0],1);
Vmath::Vvtvp (nqtot,&g5[0],1,&wsp3[0],1,&wsp8[0],1,&wsp8[0],1);
Vmath::Vvtvvtp(nqtot,&g2[0],1,&wsp3[0],1,&g4[0],1,&wsp1[0],1,&wsp9[0],1);
Vmath::Vvtvp (nqtot,&g5[0],1,&wsp2[0],1,&wsp9[0],1,&wsp9[0],1);
// Step 4.
// Multiply by quadrature metric
// Perform inner product w.r.t derivative bases.
IProductWRTBase_SumFacKernel(dbase0,base1,base2,wsp7,wsp1, wsp,false,true,true);
IProductWRTBase_SumFacKernel(base0,dbase1,base2,wsp8,wsp2, wsp,true,false,true);
IProductWRTBase_SumFacKernel(base0,base1,dbase2,wsp9,outarray,wsp,true,true,false);
// Step 5.
// Sum contributions from wsp1, wsp2 and outarray.
Vmath::Vadd(m_ncoeffs,wsp1.get(),1,outarray.get(),1,outarray.get(),1);
Vmath::Vadd(m_ncoeffs,wsp2.get(),1,outarray.get(),1,outarray.get(),1);
}
void Nektar::LocalRegions::PrismExp::v_MassMatrixOp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const StdRegions::StdMatrixKey mkey 
)
protectedvirtual

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1127 of file PrismExp.cpp.

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

{
StdExpansion::MassMatrixOp_MatFree(inarray,outarray,mkey);
}
void Nektar::LocalRegions::PrismExp::v_PhysDeriv ( const Array< OneD, const NekDouble > &  u_physical,
Array< OneD, NekDouble > &  out_dxi1,
Array< OneD, NekDouble > &  out_dxi2,
Array< OneD, NekDouble > &  out_dxi3 
)
protectedvirtual

Calculate the derivative of the physical points.

The derivative is evaluated at the nodal physical points. Derivatives with respect to the local Cartesian coordinates.

$\begin{Bmatrix} \frac {\partial} {\partial \xi_1} \\ \frac {\partial} {\partial \xi_2} \\ \frac {\partial} {\partial \xi_3} \end{Bmatrix} = \begin{Bmatrix} \frac 2 {(1-\eta_3)} \frac \partial {\partial \bar \eta_1} \\ \frac {\partial} {\partial \xi_2} \ \ \frac {(1 + \bar \eta_1)} {(1 - \eta_3)} \frac \partial {\partial \bar \eta_1} + \frac {\partial} {\partial \eta_3} \end{Bmatrix}$

Reimplemented from Nektar::StdRegions::StdPrismExp.

Definition at line 135 of file PrismExp.cpp.

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

{
int nqtot = GetTotPoints();
Array<TwoD, const NekDouble> df =
m_metricinfo->GetDerivFactors(GetPointsKeys());
Array<OneD, NekDouble> diff0(nqtot);
Array<OneD, NekDouble> diff1(nqtot);
Array<OneD, NekDouble> diff2(nqtot);
StdPrismExp::v_PhysDeriv(inarray, diff0, diff1, diff2);
{
if(out_d0.num_elements())
{
Vmath::Vmul (nqtot,&df[0][0],1,&diff0[0],1,&out_d0[0],1);
Vmath::Vvtvp (nqtot,&df[1][0],1,&diff1[0],1,&out_d0[0],1,&out_d0[0],1);
Vmath::Vvtvp (nqtot,&df[2][0],1,&diff2[0],1,&out_d0[0],1,&out_d0[0],1);
}
if(out_d1.num_elements())
{
Vmath::Vmul (nqtot,&df[3][0],1,&diff0[0],1,&out_d1[0],1);
Vmath::Vvtvp (nqtot,&df[4][0],1,&diff1[0],1,&out_d1[0],1,&out_d1[0],1);
Vmath::Vvtvp (nqtot,&df[5][0],1,&diff2[0],1,&out_d1[0],1,&out_d1[0],1);
}
if(out_d2.num_elements())
{
Vmath::Vmul (nqtot,&df[6][0],1,&diff0[0],1,&out_d2[0],1);
Vmath::Vvtvp (nqtot,&df[7][0],1,&diff1[0],1,&out_d2[0],1,&out_d2[0],1);
Vmath::Vvtvp (nqtot,&df[8][0],1,&diff2[0],1,&out_d2[0],1,&out_d2[0],1);
}
}
else // regular geometry
{
if(out_d0.num_elements())
{
Vmath::Smul (nqtot,df[0][0],&diff0[0],1,&out_d0[0],1);
Blas::Daxpy (nqtot,df[1][0],&diff1[0],1,&out_d0[0],1);
Blas::Daxpy (nqtot,df[2][0],&diff2[0],1,&out_d0[0],1);
}
if(out_d1.num_elements())
{
Vmath::Smul (nqtot,df[3][0],&diff0[0],1,&out_d1[0],1);
Blas::Daxpy (nqtot,df[4][0],&diff1[0],1,&out_d1[0],1);
Blas::Daxpy (nqtot,df[5][0],&diff2[0],1,&out_d1[0],1);
}
if(out_d2.num_elements())
{
Vmath::Smul (nqtot,df[6][0],&diff0[0],1,&out_d2[0],1);
Blas::Daxpy (nqtot,df[7][0],&diff1[0],1,&out_d2[0],1);
Blas::Daxpy (nqtot,df[8][0],&diff2[0],1,&out_d2[0],1);
}
}
}
NekDouble Nektar::LocalRegions::PrismExp::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.

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 #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::StdExpansion3D.

Definition at line 479 of file PrismExp.cpp.

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

Referenced by v_StdPhysEvaluate().

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

Reimplemented from Nektar::StdRegions::StdExpansion.

Definition at line 1183 of file PrismExp.cpp.

References Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::eOrtho_B, Nektar::StdRegions::StdExpansion::m_base, Nektar::StdRegions::StdExpansion::m_ncoeffs, Vmath::Vcopy(), and Vmath::Zero().

{
int n_coeffs = inarray.num_elements();
int nquad0 = m_base[0]->GetNumPoints();
int nquad1 = m_base[1]->GetNumPoints();
int nquad2 = m_base[2]->GetNumPoints();
int nqtot = nquad0*nquad1*nquad2;
int nmodes0 = m_base[0]->GetNumModes();
int nmodes1 = m_base[1]->GetNumModes();
int nmodes2 = m_base[2]->GetNumModes();
int numMax = nmodes0;
Array<OneD, NekDouble> coeff (n_coeffs);
Array<OneD, NekDouble> coeff_tmp1(n_coeffs, 0.0);
Array<OneD, NekDouble> coeff_tmp2(n_coeffs, 0.0);
Array<OneD, NekDouble> phys_tmp(nqtot,0.0);
Array<OneD, NekDouble> tmp, tmp2, tmp3, tmp4;
Vmath::Vcopy(n_coeffs,inarray,1,coeff_tmp2,1);
const LibUtilities::PointsKey Pkey0 = m_base[0]->GetPointsKey();
const LibUtilities::PointsKey Pkey1 = m_base[1]->GetPointsKey();
const LibUtilities::PointsKey Pkey2 = m_base[2]->GetPointsKey();
LibUtilities::BasisKey b0 = m_base[0]->GetBasisKey();
LibUtilities::BasisKey b1 = m_base[1]->GetBasisKey();
LibUtilities::BasisKey b2 = m_base[2]->GetBasisKey();
LibUtilities::BasisKey bortho0(
LibUtilities::eOrtho_A, nmodes0, Pkey0);
LibUtilities::BasisKey bortho1(
LibUtilities::eOrtho_A, nmodes1, Pkey1);
LibUtilities::BasisKey bortho2(
LibUtilities::eOrtho_B, nmodes2, Pkey2);
Vmath::Zero(n_coeffs, coeff_tmp2, 1);
int cnt = 0;
m_PrismExp = MemoryManager<StdRegions::StdPrismExp>
::AllocateSharedPtr(b0, b1, b2);
m_OrthoPrismExp = MemoryManager<StdRegions::StdPrismExp>
::AllocateSharedPtr(bortho0, bortho1, bortho2);
m_PrismExp ->BwdTrans(inarray,phys_tmp);
m_OrthoPrismExp->FwdTrans(phys_tmp, coeff);
Vmath::Zero(m_ncoeffs,outarray,1);
// filtering
for (int u = 0; u < numMax; ++u)
{
for (int i = 0; i < numMax-u; ++i)
{
Vmath::Vcopy(numMin-u,
tmp = coeff+cnt,1,
tmp2 = coeff_tmp1+cnt,1);
cnt += numMax - u;
}
}
m_OrthoPrismExp->BwdTrans(coeff_tmp1,phys_tmp);
m_PrismExp ->FwdTrans(phys_tmp, outarray);
}
NekDouble Nektar::LocalRegions::PrismExp::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 471 of file PrismExp.cpp.

References v_PhysEvaluate().

{
// Evaluate point in local (eta) coordinates.
return StdPrismExp::v_PhysEvaluate(Lcoord,physvals);
}

Member Data Documentation

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

Definition at line 207 of file PrismExp.h.

Referenced by CreateMatrix(), v_FwdTrans(), and v_GetLocMatrix().

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

Definition at line 208 of file PrismExp.h.

Referenced by v_DropLocStaticCondMatrix(), and v_GetLocStaticCondMatrix().