Nektar++
Public Member Functions | Protected Member Functions | Private Member Functions | List of all members
Nektar::MultiRegions::ExpList2D Class Reference

Abstraction of a two-dimensional multi-elemental expansion which is merely a collection of local expansions. More...

#include <ExpList2D.h>

Inheritance diagram for Nektar::MultiRegions::ExpList2D:
Inheritance graph
[legend]
Collaboration diagram for Nektar::MultiRegions::ExpList2D:
Collaboration graph
[legend]

Public Member Functions

 ExpList2D ()
 Default constructor. More...
 
 ExpList2D (const ExpList2D &In, const bool DeclareCoeffPhysArrays=true)
 Copy constructor. More...
 
 ExpList2D (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &graph2D, const bool DelcareCoeffPhysArrays=true, const std::string &var="DefaultVar")
 Sets up a list of local expansions based on an input mesh. More...
 
 ExpList2D (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::ExpansionMap &expansions, const bool DeclareCoeffPhysArrays=true)
 Sets up a list of local expansions based on an expansion Map. More...
 
 ExpList2D (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::BasisKey &TriBa, const LibUtilities::BasisKey &TriBb, const LibUtilities::BasisKey &QuadBa, const LibUtilities::BasisKey &QuadBb, const SpatialDomains::MeshGraphSharedPtr &graph2D, const LibUtilities::PointsType TriNb=LibUtilities::SIZE_PointsType)
 Sets up a list of local expansions based on an input mesh and separately defined basiskeys. More...
 
 ExpList2D (const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, const ExpListSharedPtr > &bndConstraint, const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > &bndCond, const LocalRegions::ExpansionVector &locexp, const SpatialDomains::MeshGraphSharedPtr &graph3D, const PeriodicMap &periodicFaces, const bool DeclareCoeffPhysArrays=true, const std::string variable="DefaultVar")
 
 ExpList2D (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::CompositeMap &domain, const SpatialDomains::MeshGraphSharedPtr &graph3D, const std::string variable="DefaultVar")
 Specialised constructor for Neumann boundary conditions in DisContField3D and ContField3D. More...
 
virtual ~ExpList2D ()
 Destructor. More...
 
- Public Member Functions inherited from Nektar::MultiRegions::ExpList
 ExpList ()
 The default constructor. More...
 
 ExpList (const LibUtilities::SessionReaderSharedPtr &pSession)
 The default constructor. More...
 
 ExpList (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 The default constructor. More...
 
 ExpList (const ExpList &in, const bool DeclareCoeffPhysArrays=true)
 The copy constructor. More...
 
virtual ~ExpList ()
 The default destructor. More...
 
int GetNcoeffs (void) const
 Returns the total number of local degrees of freedom $N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_m$. More...
 
int GetNcoeffs (const int eid) const
 Returns the total number of local degrees of freedom for element eid. More...
 
ExpansionType GetExpType (void)
 Returns the type of the expansion. More...
 
void SetExpType (ExpansionType Type)
 Returns the type of the expansion. More...
 
int EvalBasisNumModesMax (void) const
 Evaulates the maximum number of modes in the elemental basis order over all elements. More...
 
const Array< OneD, int > EvalBasisNumModesMaxPerExp (void) const
 Returns the vector of the number of modes in the elemental basis order over all elements. More...
 
int GetTotPoints (void) const
 Returns the total number of quadrature points m_npoints $=Q_{\mathrm{tot}}$. More...
 
int GetTotPoints (const int eid) const
 Returns the total number of quadrature points for eid's element $=Q_{\mathrm{tot}}$. More...
 
int GetNpoints (void) const
 Returns the total number of quadrature points m_npoints $=Q_{\mathrm{tot}}$. More...
 
int Get1DScaledTotPoints (const NekDouble scale) const
 Returns the total number of qudature points scaled by the factor scale on each 1D direction. More...
 
void SetWaveSpace (const bool wavespace)
 Sets the wave space to the one of the possible configuration true or false. More...
 
void SetModifiedBasis (const bool modbasis)
 Set Modified Basis for the stability analysis. More...
 
void SetPhys (int i, NekDouble val)
 Set the i th value of m_phys to value val. More...
 
bool GetWaveSpace (void) const
 This function returns the third direction expansion condition, which can be in wave space (coefficient) or not It is stored in the variable m_WaveSpace. More...
 
void SetPhys (const Array< OneD, const NekDouble > &inarray)
 Fills the array m_phys. More...
 
void SetPhysArray (Array< OneD, NekDouble > &inarray)
 Sets the array m_phys. More...
 
void SetPhysState (const bool physState)
 This function manually sets whether the array of physical values $\boldsymbol{u}_l$ (implemented as m_phys) is filled or not. More...
 
bool GetPhysState (void) const
 This function indicates whether the array of physical values $\boldsymbol{u}_l$ (implemented as m_phys) is filled or not. More...
 
NekDouble PhysIntegral (void)
 This function integrates a function $f(\boldsymbol{x})$ over the domain consisting of all the elements of the expansion. More...
 
NekDouble PhysIntegral (const Array< OneD, const NekDouble > &inarray)
 This function integrates a function $f(\boldsymbol{x})$ over the domain consisting of all the elements of the expansion. More...
 
void IProductWRTBase_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function calculates the inner product of a function $f(\boldsymbol{x})$ with respect to all {local} expansion modes $\phi_n^e(\boldsymbol{x})$. More...
 
void IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
 
void IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function calculates the inner product of a function $f(\boldsymbol{x})$ with respect to the derivative (in direction. More...
 
void IProductWRTDerivBase (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &outarray)
 This function calculates the inner product of a function $f(\boldsymbol{x})$ with respect to the derivative (in direction. More...
 
void FwdTrans_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function elementally evaluates the forward transformation of a function $u(\boldsymbol{x})$ onto the global spectral/hp expansion. More...
 
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
 
void MultiplyByElmtInvMass (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function elementally mulplies the coefficient space of Sin my the elemental inverse of the mass matrix. More...
 
void MultiplyByInvMassMatrix (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
 
void SmoothField (Array< OneD, NekDouble > &field)
 Smooth a field across elements. More...
 
void HelmSolve (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const FlagList &flags, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff=StdRegions::NullVarCoeffMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 Solve helmholtz problem. More...
 
void LinearAdvectionDiffusionReactionSolve (const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 Solve Advection Diffusion Reaction. More...
 
void LinearAdvectionReactionSolve (const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 Solve Advection Diffusion Reaction. More...
 
void FwdTrans_BndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void BwdTrans_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function elementally evaluates the backward transformation of the global spectral/hp element expansion. More...
 
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
 
void GetCoords (Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1=NullNekDouble1DArray, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
 This function calculates the coordinates of all the elemental quadrature points $\boldsymbol{x}_i$. More...
 
void HomogeneousFwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
 
void HomogeneousBwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
 
void DealiasedProd (const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
 
void GetBCValues (Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
 
void NormVectorIProductWRTBase (Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
 
void ApplyGeomInfo ()
 Apply geometry information to each expansion. More...
 
void Reset ()
 Reset geometry information and reset matrices. More...
 
void WriteTecplotHeader (std::ostream &outfile, std::string var="")
 
void WriteTecplotZone (std::ostream &outfile, int expansion=-1)
 
void WriteTecplotField (std::ostream &outfile, int expansion=-1)
 
void WriteTecplotConnectivity (std::ostream &outfile, int expansion=-1)
 
void WriteVtkHeader (std::ostream &outfile)
 
void WriteVtkFooter (std::ostream &outfile)
 
void WriteVtkPieceHeader (std::ostream &outfile, int expansion, int istrip=0)
 
void WriteVtkPieceFooter (std::ostream &outfile, int expansion)
 
void WriteVtkPieceData (std::ostream &outfile, int expansion, std::string var="v")
 
int GetCoordim (int eid)
 This function returns the dimension of the coordinates of the element eid. More...
 
void SetCoeff (int i, NekDouble val)
 Set the i th coefficiient in m_coeffs to value val. More...
 
void SetCoeffs (int i, NekDouble val)
 Set the i th coefficiient in m_coeffs to value val. More...
 
void SetCoeffsArray (Array< OneD, NekDouble > &inarray)
 Set the m_coeffs array to inarray. More...
 
const Array< OneD, const NekDouble > & GetCoeffs () const
 This function returns (a reference to) the array $\boldsymbol{\hat{u}}_l$ (implemented as m_coeffs) containing all local expansion coefficients. More...
 
void ImposeDirichletConditions (Array< OneD, NekDouble > &outarray)
 Impose Dirichlet Boundary Conditions onto Array. More...
 
void FillBndCondFromField (void)
 Fill Bnd Condition expansion from the values stored in expansion. More...
 
void LocalToGlobal (void)
 Put the coefficients into global ordering using m_coeffs. More...
 
void GlobalToLocal (void)
 Put the coefficients into local ordering and place in m_coeffs. More...
 
NekDouble GetCoeff (int i)
 Get the i th value (coefficient) of m_coeffs. More...
 
NekDouble GetCoeffs (int i)
 Get the i th value (coefficient) of m_coeffs. More...
 
const Array< OneD, const NekDouble > & GetPhys () const
 This function returns (a reference to) the array $\boldsymbol{u}_l$ (implemented as m_phys) containing the function $u^{\delta}(\boldsymbol{x})$ evaluated at the quadrature points. More...
 
NekDouble Linf (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
 This function calculates the $L_\infty$ error of the global spectral/hp element approximation. More...
 
NekDouble L2 (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
 This function calculates the $L_2$ error with respect to soln of the global spectral/hp element approximation. More...
 
NekDouble H1 (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
 Calculates the $H^1$ error of the global spectral/hp element approximation. More...
 
NekDouble Integral (const Array< OneD, const NekDouble > &inarray)
 
Array< OneD, const NekDoubleHomogeneousEnergy (void)
 This function calculates the energy associated with each one of the modesof a 3D homogeneous nD expansion. More...
 
void SetHomo1DSpecVanVisc (Array< OneD, NekDouble > visc)
 This function sets the Spectral Vanishing Viscosity in homogeneous1D expansion. More...
 
Array< OneD, const unsigned int > GetZIDs (void)
 This function returns a vector containing the wave numbers in z-direction associated with the 3D homogenous expansion. Required if a parellelisation is applied in the Fourier direction. More...
 
LibUtilities::TranspositionSharedPtr GetTransposition (void)
 This function returns the transposition class associaed with the homogeneous expansion. More...
 
NekDouble GetHomoLen (void)
 This function returns the Width of homogeneous direction associaed with the homogeneous expansion. More...
 
Array< OneD, const unsigned int > GetYIDs (void)
 This function returns a vector containing the wave numbers in y-direction associated with the 3D homogenous expansion. Required if a parellelisation is applied in the Fourier direction. More...
 
void PhysInterp1DScaled (const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function interpolates the physical space points in inarray to outarray using the same points defined in the expansion but where the number of points are rescaled by 1DScale. More...
 
void PhysGalerkinProjection1DScaled (const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function Galerkin projects the physical space points in inarray to outarray where inarray is assumed to be defined in the expansion but where the number of points are rescaled by 1DScale. More...
 
int GetExpSize (void)
 This function returns the number of elements in the expansion. More...
 
int GetNumElmts (void)
 This function returns the number of elements in the expansion which may be different for a homogeoenous extended expansionp. More...
 
const boost::shared_ptr< LocalRegions::ExpansionVectorGetExp () const
 This function returns the vector of elements in the expansion. More...
 
LocalRegions::ExpansionSharedPtrGetExp (int n) const
 This function returns (a shared pointer to) the local elemental expansion of the $n^{\mathrm{th}}$ element. More...
 
LocalRegions::ExpansionSharedPtrGetExp (const Array< OneD, const NekDouble > &gloCoord)
 This function returns (a shared pointer to) the local elemental expansion containing the arbitrary point given by gloCoord. More...
 
int GetExpIndex (const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0, bool returnNearestElmt=false)
 
int GetExpIndex (const Array< OneD, const NekDouble > &gloCoords, Array< OneD, NekDouble > &locCoords, NekDouble tol=0.0, bool returnNearestElmt=false)
 
int GetCoeff_Offset (int n) const
 Get the start offset position for a global list of m_coeffs correspoinding to element n. More...
 
int GetPhys_Offset (int n) const
 Get the start offset position for a global list of m_phys correspoinding to element n. More...
 
int GetOffset_Elmt_Id (int n) const
 Get the element id associated with the n th consecutive block of data in m_phys and m_coeffs. More...
 
Array< OneD, NekDouble > & UpdateCoeffs ()
 This function returns (a reference to) the array $\boldsymbol{\hat{u}}_l$ (implemented as m_coeffs) containing all local expansion coefficients. More...
 
Array< OneD, NekDouble > & UpdatePhys ()
 This function returns (a reference to) the array $\boldsymbol{u}_l$ (implemented as m_phys) containing the function $u^{\delta}(\boldsymbol{x})$ evaluated at the quadrature points. More...
 
void PhysDeriv (Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
 
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)
 This function discretely evaluates the derivative of a function $f(\boldsymbol{x})$ on the domain consisting of all elements of the expansion. More...
 
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
 
const Array< OneD, const boost::shared_ptr< ExpList > > & GetBndCondExpansions ()
 
boost::shared_ptr< ExpList > & UpdateBndCondExpansion (int i)
 
void Upwind (const Array< OneD, const Array< OneD, NekDouble > > &Vec, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
 
void Upwind (const Array< OneD, const NekDouble > &Vn, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
 
boost::shared_ptr< ExpList > & GetTrace ()
 
boost::shared_ptr< AssemblyMapDG > & GetTraceMap (void)
 
const Array< OneD, const int > & GetTraceBndMap (void)
 
void GetNormals (Array< OneD, Array< OneD, NekDouble > > &normals)
 
void AddTraceIntegral (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
void AddTraceIntegral (const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
void AddFwdBwdTraceIntegral (const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
 
void GetFwdBwdTracePhys (Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
 
void GetFwdBwdTracePhys (const Array< OneD, const NekDouble > &field, Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
 
void ExtractTracePhys (Array< OneD, NekDouble > &outarray)
 
void ExtractTracePhys (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & GetBndConditions ()
 
Array< OneD, SpatialDomains::BoundaryConditionShPtr > & UpdateBndConditions ()
 
void EvaluateBoundaryConditions (const NekDouble time=0.0, const std::string varName="", const NekDouble=NekConstants::kNekUnsetDouble, const NekDouble=NekConstants::kNekUnsetDouble)
 
void GeneralMatrixOp (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
 This function calculates the result of the multiplication of a matrix of type specified by mkey with a vector given by inarray. More...
 
void GeneralMatrixOp_IterPerExp (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void SetUpPhysNormals ()
 
void GetBoundaryToElmtMap (Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
 
void GeneralGetFieldDefinitions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef, int NumHomoDir=0, int NumHomoStrip=1, Array< OneD, LibUtilities::BasisSharedPtr > &HomoBasis=LibUtilities::NullBasisSharedPtr1DArray, std::vector< NekDouble > &HomoLen=LibUtilities::NullNekDoubleVector, std::vector< unsigned int > &HomoZIDs=LibUtilities::NullUnsignedIntVector, std::vector< unsigned int > &HomoYIDs=LibUtilities::NullUnsignedIntVector)
 
const NekOptimize::GlobalOptParamSharedPtrGetGlobalOptParam (void)
 
map< int, RobinBCInfoSharedPtrGetRobinBCInfo ()
 
void GetPeriodicEntities (PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces=NullPeriodicMap)
 
std::vector< LibUtilities::FieldDefinitionsSharedPtrGetFieldDefinitions ()
 
void GetFieldDefinitions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef)
 
void AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
 Append the element data listed in elements fielddef->m_ElementIDs onto fielddata. More...
 
void AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, Array< OneD, NekDouble > &coeffs)
 Append the data in coeffs listed in elements fielddef->m_ElementIDs onto fielddata. More...
 
void ExtractElmtDataToCoeffs (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
 Extract the data in fielddata into the coeffs using the basic ExpList Elemental expansions rather than planes in homogeneous case. More...
 
void ExtractCoeffsToCoeffs (const boost::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
 Extract the data from fromField using fromExpList the coeffs using the basic ExpList Elemental expansions rather than planes in homogeneous case. More...
 
void ExtractDataToCoeffs (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
 Extract the data in fielddata into the coeffs. More...
 
boost::shared_ptr< ExpListGetSharedThisPtr ()
 Returns a shared pointer to the current object. More...
 
boost::shared_ptr< LibUtilities::SessionReaderGetSession ()
 Returns the session object. More...
 
boost::shared_ptr< LibUtilities::CommGetComm ()
 Returns the comm object. More...
 
SpatialDomains::MeshGraphSharedPtr GetGraph ()
 
LibUtilities::BasisSharedPtr GetHomogeneousBasis (void)
 
boost::shared_ptr< ExpList > & GetPlane (int n)
 
void CreateCollections (Collections::ImplementationType ImpType=Collections::eNoImpType)
 Construct collections of elements containing a single element type and polynomial order from the list of expansions. More...
 

Protected Member Functions

void v_Upwind (const Array< OneD, const NekDouble > &Vn, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
 Upwind the Fwd and Bwd states based on the one- dimensional normal velocity field given by Vn. More...
 
void v_GetNormals (Array< OneD, Array< OneD, NekDouble > > &normals)
 For each local element, copy the normals stored in the element list into the array normals. More...
 
- Protected Member Functions inherited from Nektar::MultiRegions::ExpList
boost::shared_ptr< DNekMatGenGlobalMatrixFull (const GlobalLinSysKey &mkey, const boost::shared_ptr< AssemblyMapCG > &locToGloMap)
 
const DNekScalBlkMatSharedPtr GenBlockMatrix (const GlobalMatrixKey &gkey)
 This function assembles the block diagonal matrix of local matrices of the type mtype. More...
 
const DNekScalBlkMatSharedPtrGetBlockMatrix (const GlobalMatrixKey &gkey)
 
void MultiplyByBlockMatrix (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
boost::shared_ptr< GlobalMatrixGenGlobalMatrix (const GlobalMatrixKey &mkey, const boost::shared_ptr< AssemblyMapCG > &locToGloMap)
 Generates a global matrix from the given key and map. More...
 
void GlobalEigenSystem (const boost::shared_ptr< DNekMat > &Gmat, Array< OneD, NekDouble > &EigValsReal, Array< OneD, NekDouble > &EigValsImag, Array< OneD, NekDouble > &EigVecs=NullNekDouble1DArray)
 
boost::shared_ptr< GlobalLinSysGenGlobalLinSys (const GlobalLinSysKey &mkey, const boost::shared_ptr< AssemblyMapCG > &locToGloMap)
 This operation constructs the global linear system of type mkey. More...
 
boost::shared_ptr< GlobalLinSysGenGlobalBndLinSys (const GlobalLinSysKey &mkey, const AssemblyMapSharedPtr &locToGloMap)
 Generate a GlobalLinSys from information provided by the key "mkey" and the mapping provided in LocToGloBaseMap. More...
 
void ReadGlobalOptimizationParameters ()
 
virtual int v_GetNumElmts (void)
 
virtual const Array< OneD, const boost::shared_ptr< ExpList > > & v_GetBndCondExpansions (void)
 
virtual boost::shared_ptr< ExpList > & v_UpdateBndCondExpansion (int i)
 
virtual void v_Upwind (const Array< OneD, const Array< OneD, NekDouble > > &Vec, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
 
virtual boost::shared_ptr< ExpList > & v_GetTrace ()
 
virtual boost::shared_ptr< AssemblyMapDG > & v_GetTraceMap ()
 
virtual const Array< OneD, const int > & v_GetTraceBndMap ()
 
virtual void v_AddTraceIntegral (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddTraceIntegral (const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddFwdBwdTraceIntegral (const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
 
virtual void v_GetFwdBwdTracePhys (Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
 
virtual void v_GetFwdBwdTracePhys (const Array< OneD, const NekDouble > &field, Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
 
virtual void v_ExtractTracePhys (Array< OneD, NekDouble > &outarray)
 
virtual void v_ExtractTracePhys (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_MultiplyByInvMassMatrix (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
 
virtual void v_HelmSolve (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const FlagList &flags, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff, const Array< OneD, const NekDouble > &dirForcing)
 
virtual void v_LinearAdvectionDiffusionReactionSolve (const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 
virtual void v_LinearAdvectionReactionSolve (const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 
virtual void v_ImposeDirichletConditions (Array< OneD, NekDouble > &outarray)
 
virtual void v_FillBndCondFromField ()
 
virtual void v_Reset ()
 Reset geometry information, metrics, matrix managers and geometry information. More...
 
virtual void v_LocalToGlobal (void)
 
virtual void v_GlobalToLocal (void)
 
virtual void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
 
virtual void v_BwdTrans_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
 
virtual void v_FwdTrans_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_SmoothField (Array< OneD, NekDouble > &field)
 
virtual void v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
 
virtual void v_IProductWRTBase_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_GeneralMatrixOp (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
 
virtual void v_GetCoords (Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
 
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)
 
virtual void v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
 
virtual void v_PhysDeriv (Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
 
virtual void v_HomogeneousFwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
 
virtual void v_HomogeneousBwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
 
virtual void v_DealiasedProd (const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
 
virtual void v_GetBCValues (Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
 
virtual void v_NormVectorIProductWRTBase (Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
 
virtual void v_GetBoundaryToElmtMap (Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
 
virtual std::vector< LibUtilities::FieldDefinitionsSharedPtrv_GetFieldDefinitions (void)
 
virtual void v_GetFieldDefinitions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef)
 
virtual void v_AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
 
virtual void v_AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, Array< OneD, NekDouble > &coeffs)
 
virtual void v_ExtractDataToCoeffs (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
 Extract data from raw field data into expansion list. More...
 
virtual void v_ExtractCoeffsToCoeffs (const boost::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
 
virtual void v_WriteTecplotHeader (std::ostream &outfile, std::string var="")
 
virtual void v_WriteTecplotZone (std::ostream &outfile, int expansion)
 
virtual void v_WriteTecplotField (std::ostream &outfile, int expansion)
 
virtual void v_WriteTecplotConnectivity (std::ostream &outfile, int expansion)
 
virtual void v_WriteVtkPieceData (std::ostream &outfile, int expansion, std::string var)
 
virtual NekDouble v_L2 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
 
virtual NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray)
 
virtual Array< OneD, const NekDoublev_HomogeneousEnergy (void)
 
virtual LibUtilities::TranspositionSharedPtr v_GetTransposition (void)
 
virtual NekDouble v_GetHomoLen (void)
 
virtual Array< OneD, const unsigned int > v_GetZIDs (void)
 
virtual Array< OneD, const unsigned int > v_GetYIDs (void)
 
void ExtractFileBCs (const std::string &fileName, const std::string &varName, const boost::shared_ptr< ExpList > locExpList)
 

Private Member Functions

void SetCoeffPhysOffsets (void)
 Definition of the total number of degrees of freedom and quadrature points and offsets to access datax. More...
 
virtual void v_SetUpPhysNormals ()
 Set up the normals on each expansion. More...
 
virtual void v_ReadGlobalOptimizationParameters ()
 
virtual void v_WriteVtkPieceHeader (std::ostream &outfile, int expansion, int istrip)
 
virtual void v_PhysInterp1DScaled (const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_PhysGalerkinProjection1DScaled (const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 

Additional Inherited Members

- Public Attributes inherited from Nektar::MultiRegions::ExpList
ExpansionType m_expType
 
- Static Protected Member Functions inherited from Nektar::MultiRegions::ExpList
static SpatialDomains::BoundaryConditionShPtr GetBoundaryCondition (const SpatialDomains::BoundaryConditionCollection &collection, unsigned int index, const std::string &variable)
 
- Protected Attributes inherited from Nektar::MultiRegions::ExpList
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
LibUtilities::SessionReaderSharedPtr m_session
 Session. More...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 Mesh associated with this expansion list. More...
 
int m_ncoeffs
 The total number of local degrees of freedom. m_ncoeffs $=N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_l$. More...
 
int m_npoints
 
Array< OneD, NekDoublem_coeffs
 Concatenation of all local expansion coefficients. More...
 
Array< OneD, NekDoublem_phys
 The global expansion evaluated at the quadrature points. More...
 
bool m_physState
 The state of the array m_phys. More...
 
boost::shared_ptr< LocalRegions::ExpansionVectorm_exp
 The list of local expansions. More...
 
Collections::CollectionVector m_collections
 
std::vector< int > m_coll_coeff_offset
 Offset of elemental data into the array m_coeffs. More...
 
std::vector< int > m_coll_phys_offset
 Offset of elemental data into the array m_phys. More...
 
Array< OneD, int > m_coeff_offset
 Offset of elemental data into the array m_coeffs. More...
 
Array< OneD, int > m_phys_offset
 Offset of elemental data into the array m_phys. More...
 
Array< OneD, int > m_offset_elmt_id
 Array containing the element id m_offset_elmt_id[n] that the n^th consecutive block of data in m_coeffs and m_phys is associated, i.e. for an array of constant expansion size and single shape elements m_phys[n*m_npoints] is the data related to m_exp[m_offset_elmt_id[n]];. More...
 
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
 
BlockMatrixMapShPtr m_blockMat
 
bool m_WaveSpace
 

Detailed Description

Abstraction of a two-dimensional multi-elemental expansion which is merely a collection of local expansions.

This multi-elemental expansion, which does not exhibit any coupling between the expansion on the separate elements, can be formulated as,

\[u^{\delta}(\boldsymbol{x}_i)=\sum_{e=1}^{{N_{\mathrm{el}}}} \sum_{n=0}^{N^{e}_m-1}\hat{u}_n^e\phi_n^e(\boldsymbol{x}_i).\]

where ${N_{\mathrm{el}}}$ is the number of elements and $N^{e}_m$ is the local elemental number of expansion modes. This class inherits all its variables and member functions from the base class ExpList.

Definition at line 60 of file ExpList2D.h.

Constructor & Destructor Documentation

Nektar::MultiRegions::ExpList2D::ExpList2D ( )

Default constructor.

Definition at line 68 of file ExpList2D.cpp.

References Nektar::MultiRegions::e2D, and Nektar::MultiRegions::ExpList::SetExpType().

68  :
69  ExpList()
70  {
71  SetExpType(e2D);
72  }
ExpList()
The default constructor.
Definition: ExpList.cpp:93
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:211
Nektar::MultiRegions::ExpList2D::ExpList2D ( const ExpList2D In,
const bool  DeclareCoeffPhysArrays = true 
)

Copy constructor.

Parameters
InExpList2D object to copy.

Definition at line 86 of file ExpList2D.cpp.

References Nektar::MultiRegions::e2D, and Nektar::MultiRegions::ExpList::SetExpType().

88  :
89  ExpList(In,DeclareCoeffPhysArrays)
90  {
91  SetExpType(e2D);
92  }
ExpList()
The default constructor.
Definition: ExpList.cpp:93
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:211
Nektar::MultiRegions::ExpList2D::ExpList2D ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr graph2D,
const bool  DeclareCoeffPhysArrays = true,
const std::string &  var = "DefaultVar" 
)

Sets up a list of local expansions based on an input mesh.

Given a mesh graph2D, containing information about the domain and the spectral/hp element expansion, this constructor fills the list of local expansions {m_exp} with the proper expansions, calculates the total number of quadrature points $\boldsymbol{x}_i$ and local expansion coefficients $\hat{u}^e_n$ and allocates memory for the arrays m_coeffs and m_phys.

Parameters
graph2DA mesh, containing information about the domain and the spectral/hp element expansion.

Definition at line 107 of file ExpList2D.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::MultiRegions::ExpList::CreateCollections(), Nektar::MultiRegions::e2D, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eNodalTriElec, Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::BasisKey::GetBasisType(), Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::LibUtilities::BasisKey::GetNumModes(), Nektar::LibUtilities::BasisKey::GetNumPoints(), Nektar::LibUtilities::BasisKey::GetPointsKey(), Nektar::MultiRegions::ExpList::m_coeffs, Nektar::MultiRegions::ExpList::m_globalOptParam, Nektar::MultiRegions::ExpList::m_ncoeffs, Nektar::MultiRegions::ExpList::m_npoints, Nektar::MultiRegions::ExpList::m_phys, Nektar::MultiRegions::ExpList::ReadGlobalOptimizationParameters(), SetCoeffPhysOffsets(), and Nektar::MultiRegions::ExpList::SetExpType().

111  :
112  ExpList(pSession,graph2D)
113  {
114  SetExpType(e2D);
115 
116  int elmtid=0;
122 
123  const SpatialDomains::ExpansionMap &expansions
124  = graph2D->GetExpansions(var);
125 
126  SpatialDomains::ExpansionMap::const_iterator expIt;
127  for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
128  {
130  SpatialDomains::QuadGeomSharedPtr QuadrilateralGeom;
131 
132  if ((TriangleGeom = boost::dynamic_pointer_cast<SpatialDomains
133  ::TriGeom>(expIt->second->m_geomShPtr)))
134  {
136  = expIt->second->m_basisKeyVector[0];
138  = expIt->second->m_basisKeyVector[1];
139 
140  // This is not elegantly implemented needs re-thinking.
142  {
144  TriBa.GetNumModes(),
145  TriBa.GetPointsKey());
146 
149  ::AllocateSharedPtr(newBa,TriBb,TriNb,
150  TriangleGeom);
151  Ntri->SetElmtId(elmtid++);
152  (*m_exp).push_back(Ntri);
153  }
154  else
155  {
157  ::AllocateSharedPtr(TriBa,TriBb,
158  TriangleGeom);
159  tri->SetElmtId(elmtid++);
160  (*m_exp).push_back(tri);
161  }
162  m_ncoeffs += (TriBa.GetNumModes()*(TriBa.GetNumModes()+1))/2
163  + TriBa.GetNumModes()*(TriBb.GetNumModes()
164  -TriBa.GetNumModes());
165  m_npoints += TriBa.GetNumPoints()*TriBb.GetNumPoints();
166  }
167  else if ((QuadrilateralGeom = boost::dynamic_pointer_cast<
168  SpatialDomains::QuadGeom>(expIt->second->m_geomShPtr)))
169  {
171  = expIt->second->m_basisKeyVector[0];
173  = expIt->second->m_basisKeyVector[1];
174 
176  ::AllocateSharedPtr(QuadBa,QuadBb,
177  QuadrilateralGeom);
178  quad->SetElmtId(elmtid++);
179  (*m_exp).push_back(quad);
180 
181  m_ncoeffs += QuadBa.GetNumModes()*QuadBb.GetNumModes();
182  m_npoints += QuadBa.GetNumPoints()*QuadBb.GetNumPoints();
183  }
184  else
185  {
186  ASSERTL0(false, "dynamic cast to a proper Geometry2D "
187  "failed");
188  }
189 
190  }
191 
192  // set up element numbering
193  for(int i = 0; i < (*m_exp).size(); ++i)
194  {
195  (*m_exp)[i]->SetElmtId(i);
196  }
197 
198  // Setup Default optimisation information.
199  int nel = GetExpSize();
202 
203 
204  // set up offset arrays.
206 
207  if (DeclareCoeffPhysArrays)
208  {
209  // Set up m_coeffs, m_phys.
212  }
213 
216  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
ExpList()
The default constructor.
Definition: ExpList.cpp:93
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:971
BasisType GetBasisType() const
Return type of expansion basis.
Definition: Basis.h:139
boost::shared_ptr< QuadGeom > QuadGeomSharedPtr
Definition: QuadGeom.h:54
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:926
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:909
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1837
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:887
int GetNumPoints() const
Return points order at which basis is defined.
Definition: Basis.h:128
Principle Orthogonal Functions .
Definition: BasisType.h:46
boost::shared_ptr< QuadExp > QuadExpSharedPtr
Definition: QuadExp.h:284
void SetCoeffPhysOffsets(void)
Definition of the total number of degrees of freedom and quadrature points and offsets to access data...
Definition: ExpList2D.cpp:1096
boost::shared_ptr< GeometryVector > Composite
Definition: MeshGraph.h:111
PointsKey GetPointsKey() const
Return distribution of points.
Definition: Basis.h:145
boost::shared_ptr< TriGeom > TriGeomSharedPtr
Definition: TriGeom.h:58
void CreateCollections(Collections::ImplementationType ImpType=Collections::eNoImpType)
Construct collections of elements containing a single element type and polynomial order from the list...
Definition: ExpList.cpp:2723
int GetNumModes() const
Returns the order of the basis.
Definition: Basis.h:84
Lagrange for SEM basis .
Definition: BasisType.h:53
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:211
Describes the specification for a Basis.
Definition: Basis.h:50
2D Nodal Electrostatic Points on a Triangle
Definition: PointsType.h:68
std::map< int, ExpansionShPtr > ExpansionMap
Definition: MeshGraph.h:171
boost::shared_ptr< TriExp > TriExpSharedPtr
Definition: TriExp.h:283
boost::shared_ptr< NodalTriExp > NodalTriExpSharedPtr
Definition: NodalTriExp.h:379
Nektar::MultiRegions::ExpList2D::ExpList2D ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::ExpansionMap expansions,
const bool  DeclareCoeffPhysArrays = true 
)

Sets up a list of local expansions based on an expansion Map.

Given an expansion vector expansions, containing information about the domain and the spectral/hp element expansion, this constructor fills the list of local expansions {m_exp} with the proper expansions, calculates the total number of quadrature points $\boldsymbol{x}_i$ and local expansion coefficients $\hat{u}^e_n$ and allocates memory for the arrays m_coeffs and m_phys.

Parameters
expansionsA vector containing information about the domain and the spectral/hp element expansion.

Definition at line 233 of file ExpList2D.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::MultiRegions::ExpList::CreateCollections(), Nektar::MultiRegions::e2D, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eNodalTriElec, Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::BasisKey::GetBasisType(), Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::LibUtilities::BasisKey::GetNumModes(), Nektar::LibUtilities::BasisKey::GetNumPoints(), Nektar::LibUtilities::BasisKey::GetPointsKey(), Nektar::MultiRegions::ExpList::m_coeffs, Nektar::MultiRegions::ExpList::m_globalOptParam, Nektar::MultiRegions::ExpList::m_ncoeffs, Nektar::MultiRegions::ExpList::m_npoints, Nektar::MultiRegions::ExpList::m_phys, Nektar::MultiRegions::ExpList::ReadGlobalOptimizationParameters(), SetCoeffPhysOffsets(), and Nektar::MultiRegions::ExpList::SetExpType().

236  :ExpList(pSession)
237  {
238  SetExpType(e2D);
239 
240  int elmtid=0;
246 
248  for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
249  {
251  SpatialDomains::QuadGeomSharedPtr QuadrilateralGeom;
252 
253  if ((TriangleGeom = boost::dynamic_pointer_cast<SpatialDomains
254  ::TriGeom>(expIt->second->m_geomShPtr)))
255  {
257  = expIt->second->m_basisKeyVector[0];
259  = expIt->second->m_basisKeyVector[1];
260 
261  // This is not elegantly implemented needs re-thinking.
263  {
265  TriBa.GetNumModes(),
266  TriBa.GetPointsKey());
267 
270  ::AllocateSharedPtr(newBa,TriBb,TriNb,
271  TriangleGeom);
272  Ntri->SetElmtId(elmtid++);
273  (*m_exp).push_back(Ntri);
274  }
275  else
276  {
278  ::AllocateSharedPtr(TriBa,TriBb,
279  TriangleGeom);
280  tri->SetElmtId(elmtid++);
281  (*m_exp).push_back(tri);
282  }
283  m_ncoeffs += (TriBa.GetNumModes()*(TriBa.GetNumModes()+1))/2
284  + TriBa.GetNumModes()*(TriBb.GetNumModes()
285  -TriBa.GetNumModes());
286  m_npoints += TriBa.GetNumPoints()*TriBb.GetNumPoints();
287  }
288  else if ((QuadrilateralGeom = boost::dynamic_pointer_cast<
289  SpatialDomains::QuadGeom>(expIt->second->m_geomShPtr)))
290  {
292  = expIt->second->m_basisKeyVector[0];
294  = expIt->second->m_basisKeyVector[1];
295 
297  ::AllocateSharedPtr(QuadBa,QuadBb,
298  QuadrilateralGeom);
299  quad->SetElmtId(elmtid++);
300  (*m_exp).push_back(quad);
301 
302  m_ncoeffs += QuadBa.GetNumModes()*QuadBb.GetNumModes();
303  m_npoints += QuadBa.GetNumPoints()*QuadBb.GetNumPoints();
304  }
305  else
306  {
307  ASSERTL0(false, "dynamic cast to a proper Geometry2D "
308  "failed");
309  }
310 
311  }
312 
313  // Setup Default optimisation information.
314  int nel = GetExpSize();
317 
318 
319  // set up offset arrays.
321 
322  if (DeclareCoeffPhysArrays)
323  {
324  // Set up m_coeffs, m_phys.
327  }
328 
331  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
std::map< int, ExpansionShPtr >::const_iterator ExpansionMapConstIter
Definition: MeshGraph.h:173
ExpList()
The default constructor.
Definition: ExpList.cpp:93
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:971
BasisType GetBasisType() const
Return type of expansion basis.
Definition: Basis.h:139
boost::shared_ptr< QuadGeom > QuadGeomSharedPtr
Definition: QuadGeom.h:54
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:926
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:909
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1837
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:887
int GetNumPoints() const
Return points order at which basis is defined.
Definition: Basis.h:128
Principle Orthogonal Functions .
Definition: BasisType.h:46
boost::shared_ptr< QuadExp > QuadExpSharedPtr
Definition: QuadExp.h:284
void SetCoeffPhysOffsets(void)
Definition of the total number of degrees of freedom and quadrature points and offsets to access data...
Definition: ExpList2D.cpp:1096
boost::shared_ptr< GeometryVector > Composite
Definition: MeshGraph.h:111
PointsKey GetPointsKey() const
Return distribution of points.
Definition: Basis.h:145
boost::shared_ptr< TriGeom > TriGeomSharedPtr
Definition: TriGeom.h:58
void CreateCollections(Collections::ImplementationType ImpType=Collections::eNoImpType)
Construct collections of elements containing a single element type and polynomial order from the list...
Definition: ExpList.cpp:2723
int GetNumModes() const
Returns the order of the basis.
Definition: Basis.h:84
Lagrange for SEM basis .
Definition: BasisType.h:53
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:211
Describes the specification for a Basis.
Definition: Basis.h:50
2D Nodal Electrostatic Points on a Triangle
Definition: PointsType.h:68
boost::shared_ptr< TriExp > TriExpSharedPtr
Definition: TriExp.h:283
boost::shared_ptr< NodalTriExp > NodalTriExpSharedPtr
Definition: NodalTriExp.h:379
Nektar::MultiRegions::ExpList2D::ExpList2D ( const LibUtilities::SessionReaderSharedPtr pSession,
const LibUtilities::BasisKey TriBa,
const LibUtilities::BasisKey TriBb,
const LibUtilities::BasisKey QuadBa,
const LibUtilities::BasisKey QuadBb,
const SpatialDomains::MeshGraphSharedPtr graph2D,
const LibUtilities::PointsType  TriNb = LibUtilities::SIZE_PointsType 
)

Sets up a list of local expansions based on an input mesh and separately defined basiskeys.

Given a mesh graph2D, containing information about the domain and the a list of basiskeys, this constructor fills the list of local expansions {m_exp} with the proper expansions, calculates the total number of quadrature points $\boldsymbol{x}_i$ and local expansion coefficients $\hat{u}^e_n$ and allocates memory for the arrays m_coeffs and m_phys.

Parameters
TriBaA BasisKey, containing the definition of the basis in the first coordinate direction for triangular elements
TriBbA BasisKey, containing the definition of the basis in the second coordinate direction for triangular elements
QuadBaA BasisKey, containing the definition of the basis in the first coordinate direction for quadrilateral elements
QuadBbA BasisKey, containing the definition of the basis in the second coordinate direction for quadrilateral elements
graph2DA mesh, containing information about the domain and the spectral/hp element expansion.
TriNbThe PointsType of possible nodal points

Definition at line 359 of file ExpList2D.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::MultiRegions::ExpList::CreateCollections(), Nektar::MultiRegions::e2D, Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::LibUtilities::BasisKey::GetNumModes(), Nektar::LibUtilities::BasisKey::GetNumPoints(), Nektar::MultiRegions::ExpList::m_coeffs, Nektar::MultiRegions::ExpList::m_globalOptParam, Nektar::MultiRegions::ExpList::m_ncoeffs, Nektar::MultiRegions::ExpList::m_npoints, Nektar::MultiRegions::ExpList::m_phys, Nektar::MultiRegions::ExpList::m_physState, Nektar::MultiRegions::ExpList::ReadGlobalOptimizationParameters(), SetCoeffPhysOffsets(), Nektar::MultiRegions::ExpList::SetExpType(), and Nektar::LibUtilities::SIZE_PointsType.

366  :ExpList(pSession, graph2D)
367  {
368  SetExpType(e2D);
369 
370  int elmtid=0;
375 
376  const SpatialDomains::ExpansionMap &expansions =
377  graph2D->GetExpansions();
378  m_ncoeffs = 0;
379  m_npoints = 0;
380 
381  m_physState = false;
382 
383  SpatialDomains::ExpansionMap::const_iterator expIt;
384  for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
385  {
387  SpatialDomains::QuadGeomSharedPtr QuadrilateralGeom;
388 
389  if ((TriangleGeom = boost::dynamic_pointer_cast<SpatialDomains::
390  TriGeom>(expIt->second->m_geomShPtr)))
391  {
392  if (TriNb < LibUtilities::SIZE_PointsType)
393  {
395  AllocateSharedPtr(TriBa, TriBb, TriNb,
396  TriangleGeom);
397  Ntri->SetElmtId(elmtid++);
398  (*m_exp).push_back(Ntri);
399  }
400  else
401  {
403  AllocateSharedPtr(TriBa, TriBb, TriangleGeom);
404  tri->SetElmtId(elmtid++);
405  (*m_exp).push_back(tri);
406  }
407 
408  m_ncoeffs += (TriBa.GetNumModes()*(TriBa.GetNumModes()+1))/2
409  + TriBa.GetNumModes() * (TriBb.GetNumModes() -
410  TriBa.GetNumModes());
411  m_npoints += TriBa.GetNumPoints()*TriBb.GetNumPoints();
412  }
413  else if ((QuadrilateralGeom = boost::dynamic_pointer_cast<
414  SpatialDomains::QuadGeom>(expIt->second->m_geomShPtr)))
415  {
417  AllocateSharedPtr(QuadBa, QuadBb, QuadrilateralGeom);
418  quad->SetElmtId(elmtid++);
419  (*m_exp).push_back(quad);
420 
421  m_ncoeffs += QuadBa.GetNumModes()*QuadBb.GetNumModes();
422  m_npoints += QuadBa.GetNumPoints()*QuadBb.GetNumPoints();
423  }
424  else
425  {
426  ASSERTL0(false,
427  "dynamic cast to a proper Geometry2D failed");
428  }
429 
430  }
431 
432  // Setup Default optimisation information.
433  int nel = GetExpSize();
436 
437  // Set up m_coeffs, m_phys and offset arrays.
441 
444  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
ExpList()
The default constructor.
Definition: ExpList.cpp:93
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:971
boost::shared_ptr< QuadGeom > QuadGeomSharedPtr
Definition: QuadGeom.h:54
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:926
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:909
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1837
bool m_physState
The state of the array m_phys.
Definition: ExpList.h:935
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:887
int GetNumPoints() const
Return points order at which basis is defined.
Definition: Basis.h:128
boost::shared_ptr< QuadExp > QuadExpSharedPtr
Definition: QuadExp.h:284
void SetCoeffPhysOffsets(void)
Definition of the total number of degrees of freedom and quadrature points and offsets to access data...
Definition: ExpList2D.cpp:1096
boost::shared_ptr< GeometryVector > Composite
Definition: MeshGraph.h:111
boost::shared_ptr< TriGeom > TriGeomSharedPtr
Definition: TriGeom.h:58
void CreateCollections(Collections::ImplementationType ImpType=Collections::eNoImpType)
Construct collections of elements containing a single element type and polynomial order from the list...
Definition: ExpList.cpp:2723
int GetNumModes() const
Returns the order of the basis.
Definition: Basis.h:84
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:211
std::map< int, ExpansionShPtr > ExpansionMap
Definition: MeshGraph.h:171
boost::shared_ptr< TriExp > TriExpSharedPtr
Definition: TriExp.h:283
boost::shared_ptr< NodalTriExp > NodalTriExpSharedPtr
Definition: NodalTriExp.h:379
Nektar::MultiRegions::ExpList2D::ExpList2D ( const LibUtilities::SessionReaderSharedPtr pSession,
const Array< OneD, const ExpListSharedPtr > &  bndConstraint,
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > &  bndCond,
const LocalRegions::ExpansionVector locexp,
const SpatialDomains::MeshGraphSharedPtr graph3D,
const PeriodicMap periodicFaces,
const bool  DeclareCoeffPhysArrays = true,
const std::string  variable = "DefaultVar" 
)

Specialized constructor for trace expansions. Store expansions for the trace space used in DisContField3D

Parameters
bndConstraintArray of ExpList2D objects each containing a 2D spectral/hp element expansion on a single boundary region.
bndCondArray of BoundaryCondition objects which contain information about the boundary conditions on the different boundary regions.
locexpComplete domain expansion list.
graph3D3D mesh corresponding to the expansion list.
periodicFacesList of periodic faces.
DeclareCoeffPhysArraysIf true, set up m_coeffs, m_phys arrays

Definition at line 462 of file ExpList2D.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::MultiRegions::ExpList::CreateCollections(), Nektar::MultiRegions::e2D, Nektar::SpatialDomains::eDirichlet, Nektar::LibUtilities::BasisKey::GetBasisType(), Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::LocalRegions::Expansion2D::GetGeom2D(), Nektar::LocalRegions::Expansion3D::GetGeom3D(), Nektar::StdRegions::StdExpansion::GetNfaces(), Nektar::LibUtilities::BasisKey::GetNumModes(), Nektar::LibUtilities::BasisKey::GetNumPoints(), Nektar::LibUtilities::BasisKey::GetPointsType(), Nektar::iterator, Nektar::MultiRegions::ExpList::m_coeffs, Nektar::MultiRegions::ExpList::m_globalOptParam, Nektar::MultiRegions::ExpList::m_ncoeffs, Nektar::MultiRegions::ExpList::m_npoints, Nektar::MultiRegions::ExpList::m_phys, Nektar::LibUtilities::ReduceSum, SetCoeffPhysOffsets(), Nektar::MultiRegions::ExpList::SetExpType(), and Vmath::Vsum().

470  :
471  ExpList(pSession, graph3D)
472  {
473  SetExpType(e2D);
474 
475  int i, j, id, elmtid=0;
476  set<int> facesDone;
477 
481  LocalRegions::QuadExpSharedPtr FaceQuadExp;
485 
486  // First loop over boundary conditions to renumber
487  // Dirichlet boundaries
488  for (i = 0; i < bndCond.num_elements(); ++i)
489  {
490  if (bndCond[i]->GetBoundaryConditionType()
492  {
493  for (j = 0; j < bndConstraint[i]->GetExpSize(); ++j)
494  {
495  LibUtilities::BasisKey bkey0 = bndConstraint[i]
496  ->GetExp(j)->GetBasis(0)->GetBasisKey();
497  LibUtilities::BasisKey bkey1 = bndConstraint[i]
498  ->GetExp(j)->GetBasis(1)->GetBasisKey();
499  exp2D = bndConstraint[i]->GetExp(j)
501  FaceGeom = exp2D->GetGeom2D();
502 
503  //if face is a quad
504  if((FaceQuadGeom = boost::dynamic_pointer_cast<
505  SpatialDomains::QuadGeom>(FaceGeom)))
506  {
508  ::AllocateSharedPtr(bkey0, bkey1, FaceQuadGeom);
509  facesDone.insert(FaceQuadGeom->GetFid());
510  FaceQuadExp->SetElmtId(elmtid++);
511  (*m_exp).push_back(FaceQuadExp);
512  }
513  //if face is a triangle
514  else if((FaceTriGeom = boost::dynamic_pointer_cast<
515  SpatialDomains::TriGeom>(FaceGeom)))
516  {
518  ::AllocateSharedPtr(bkey0, bkey1, FaceTriGeom);
519  facesDone.insert(FaceTriGeom->GetFid());
520  FaceTriExp->SetElmtId(elmtid++);
521  (*m_exp).push_back(FaceTriExp);
522  }
523  else
524  {
525  ASSERTL0(false,"dynamic cast to a proper face geometry failed");
526  }
527  }
528  }
529  }
530 
533  LibUtilities::BasisKey> > > faceOrders;
535  pair<LibUtilities::BasisKey,
536  LibUtilities::BasisKey> > >::iterator it;
537 
538  for(i = 0; i < locexp.size(); ++i)
539  {
540  exp3D = locexp[i]->as<LocalRegions::Expansion3D>();
541  for (j = 0; j < exp3D->GetNfaces(); ++j)
542  {
543  FaceGeom = exp3D->GetGeom3D()->GetFace(j);
544  id = FaceGeom->GetFid();
545 
546  if(facesDone.count(id) != 0)
547  {
548  continue;
549  }
550  it = faceOrders.find(id);
551 
552  if (it == faceOrders.end())
553  {
554  LibUtilities::BasisKey face_dir0
555  = locexp[i]->DetFaceBasisKey(j,0);
556  LibUtilities::BasisKey face_dir1
557  = locexp[i]->DetFaceBasisKey(j,1);
558 
559  faceOrders.insert(
560  std::make_pair(
561  id, std::make_pair(
562  FaceGeom,
563  std::make_pair(face_dir0, face_dir1))));
564  }
565  else // variable modes/points
566  {
567  LibUtilities::BasisKey face0 =
568  locexp[i]->DetFaceBasisKey(j,0);
569  LibUtilities::BasisKey face1 =
570  locexp[i]->DetFaceBasisKey(j,1);
571  LibUtilities::BasisKey existing0 =
572  it->second.second.first;
573  LibUtilities::BasisKey existing1 =
574  it->second.second.second;
575 
576  int np11 = face0 .GetNumPoints();
577  int np12 = face1 .GetNumPoints();
578  int np21 = existing0.GetNumPoints();
579  int np22 = existing1.GetNumPoints();
580  int nm11 = face0 .GetNumModes ();
581  int nm12 = face1 .GetNumModes ();
582  int nm21 = existing0.GetNumModes ();
583  int nm22 = existing1.GetNumModes ();
584 
585  if ((np22 >= np12 || np21 >= np11) &&
586  (nm22 >= nm12 || nm21 >= nm11))
587  {
588  continue;
589  }
590  else if((np22 < np12 || np21 < np11) &&
591  (nm22 < nm12 || nm21 < nm11))
592  {
593  it->second.second.first = face0;
594  it->second.second.second = face1;
595  }
596  else
597  {
598  ASSERTL0(false,
599  "inappropriate number of points/modes (max "
600  "num of points is not set with max order)");
601  }
602  }
603  }
604  }
605 
606  LibUtilities::CommSharedPtr vComm = pSession->GetComm();
607  int nproc = vComm->GetSize(); // number of processors
608  int facepr = vComm->GetRank(); // ID processor
609 
610  if (nproc > 1)
611  {
612  int fCnt = 0;
613 
614  // Count the number of faces on each partition
615  for(i = 0; i < locexp.size(); ++i)
616  {
617  fCnt += locexp[i]->GetNfaces();
618  }
619 
620  // Set up the offset and the array that will contain the list of
621  // face IDs, then reduce this across processors.
622  Array<OneD, int> faceCnt(nproc,0);
623  faceCnt[facepr] = fCnt;
624  vComm->AllReduce(faceCnt, LibUtilities::ReduceSum);
625 
626  int totFaceCnt = Vmath::Vsum(nproc, faceCnt, 1);
627  Array<OneD, int> fTotOffsets(nproc,0);
628 
629  for (i = 1; i < nproc; ++i)
630  {
631  fTotOffsets[i] = fTotOffsets[i-1] + faceCnt[i-1];
632  }
633 
634  // Local list of the edges per element
635 
636  Array<OneD, int> FacesTotID (totFaceCnt, 0);
637  Array<OneD, int> FacesTotNm0 (totFaceCnt, 0);
638  Array<OneD, int> FacesTotNm1 (totFaceCnt, 0);
639  Array<OneD, int> FacesTotPnts0(totFaceCnt, 0);
640  Array<OneD, int> FacesTotPnts1(totFaceCnt, 0);
641 
642  int cntr = fTotOffsets[facepr];
643 
644  for(i = 0; i < locexp.size(); ++i)
645  {
646  exp3D = locexp[i]->as<LocalRegions::Expansion3D>();
647 
648  int nfaces = locexp[i]->GetNfaces();
649 
650  for(j = 0; j < nfaces; ++j, ++cntr)
651  {
652  LibUtilities::BasisKey face_dir0
653  = locexp[i]->DetFaceBasisKey(j,0);
654  LibUtilities::BasisKey face_dir1
655  = locexp[i]->DetFaceBasisKey(j,1);
656 
657  FacesTotID[cntr] = exp3D->GetGeom3D()->GetFid(j);
658  FacesTotNm0[cntr] = face_dir0.GetNumModes ();
659  FacesTotNm1[cntr] = face_dir1.GetNumModes ();
660  FacesTotPnts0[cntr] = face_dir0.GetNumPoints();
661  FacesTotPnts1[cntr] = face_dir1.GetNumPoints();
662  }
663  }
664 
665  vComm->AllReduce(FacesTotID, LibUtilities::ReduceSum);
666  vComm->AllReduce(FacesTotNm0, LibUtilities::ReduceSum);
667  vComm->AllReduce(FacesTotNm1, LibUtilities::ReduceSum);
668  vComm->AllReduce(FacesTotPnts0, LibUtilities::ReduceSum);
669  vComm->AllReduce(FacesTotPnts1, LibUtilities::ReduceSum);
670 
671  for (i = 0; i < totFaceCnt; ++i)
672  {
673  it = faceOrders.find(FacesTotID[i]);
674 
675  if (it == faceOrders.end())
676  {
677  continue;
678  }
679 
680  LibUtilities::BasisKey existing0 =
681  it->second.second.first;
682  LibUtilities::BasisKey existing1 =
683  it->second.second.second;
684  LibUtilities::BasisKey face0(
685  existing0.GetBasisType(), FacesTotNm0[i],
686  LibUtilities::PointsKey(FacesTotPnts0[i],
687  existing0.GetPointsType()));
688  LibUtilities::BasisKey face1(
689  existing1.GetBasisType(), FacesTotNm1[i],
690  LibUtilities::PointsKey(FacesTotPnts1[i],
691  existing1.GetPointsType()));
692 
693  int np11 = face0 .GetNumPoints();
694  int np12 = face1 .GetNumPoints();
695  int np21 = existing0.GetNumPoints();
696  int np22 = existing1.GetNumPoints();
697  int nm11 = face0 .GetNumModes ();
698  int nm12 = face1 .GetNumModes ();
699  int nm21 = existing0.GetNumModes ();
700  int nm22 = existing1.GetNumModes ();
701 
702  if ((np22 >= np12 || np21 >= np11) &&
703  (nm22 >= nm12 || nm21 >= nm11))
704  {
705  continue;
706  }
707  else if((np22 < np12 || np21 < np11) &&
708  (nm22 < nm12 || nm21 < nm11))
709  {
710  it->second.second.first = face0;
711  it->second.second.second = face1;
712  }
713  else
714  {
715  ASSERTL0(false,
716  "inappropriate number of points/modes (max "
717  "num of points is not set with max order)");
718  }
719  }
720  }
721 
722  for (it = faceOrders.begin(); it != faceOrders.end(); ++it)
723  {
724  FaceGeom = it->second.first;
725 
726  if ((FaceQuadGeom = boost::dynamic_pointer_cast<
727  SpatialDomains::QuadGeom>(FaceGeom)))
728  {
730  ::AllocateSharedPtr(it->second.second.first,
731  it->second.second.second,
732  FaceQuadGeom);
733  FaceQuadExp->SetElmtId(elmtid++);
734  (*m_exp).push_back(FaceQuadExp);
735  }
736  else if ((FaceTriGeom = boost::dynamic_pointer_cast<
737  SpatialDomains::TriGeom>(FaceGeom)))
738  {
740  ::AllocateSharedPtr(it->second.second.first,
741  it->second.second.second,
742  FaceTriGeom);
743  FaceTriExp->SetElmtId(elmtid++);
744  (*m_exp).push_back(FaceTriExp);
745  }
746  }
747 
748  // Setup Default optimisation information.
749  int nel = GetExpSize();
750 
753 
754  // Set up offset information and array sizes
756 
757  // Set up m_coeffs, m_phys.
758  if(DeclareCoeffPhysArrays)
759  {
762  }
763 
765  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
ExpList()
The default constructor.
Definition: ExpList.cpp:93
int GetNfaces() const
This function returns the number of faces of the expansion domain.
Definition: StdExpansion.h:422
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:971
BasisType GetBasisType() const
Return type of expansion basis.
Definition: Basis.h:139
boost::shared_ptr< QuadGeom > QuadGeomSharedPtr
Definition: QuadGeom.h:54
PointsType GetPointsType() const
Return type of quadrature.
Definition: Basis.h:151
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:926
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:909
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1837
boost::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:53
SpatialDomains::Geometry2DSharedPtr GetGeom2D() const
Definition: Expansion2D.h:231
boost::shared_ptr< Expansion3D > Expansion3DSharedPtr
Definition: Expansion2D.h:48
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:887
int GetNumPoints() const
Return points order at which basis is defined.
Definition: Basis.h:128
Defines a specification for a set of points.
Definition: Points.h:58
boost::shared_ptr< QuadExp > QuadExpSharedPtr
Definition: QuadExp.h:284
void SetCoeffPhysOffsets(void)
Definition of the total number of degrees of freedom and quadrature points and offsets to access data...
Definition: ExpList2D.cpp:1096
boost::shared_ptr< Geometry2D > Geometry2DSharedPtr
Definition: Geometry2D.h:59
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
boost::shared_ptr< TriGeom > TriGeomSharedPtr
Definition: TriGeom.h:58
void CreateCollections(Collections::ImplementationType ImpType=Collections::eNoImpType)
Construct collections of elements containing a single element type and polynomial order from the list...
Definition: ExpList.cpp:2723
const LibUtilities::BasisKey DetFaceBasisKey(const int i, const int k) const
Definition: StdExpansion.h:324
int GetNumModes() const
Returns the order of the basis.
Definition: Basis.h:84
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:211
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:714
Describes the specification for a Basis.
Definition: Basis.h:50
boost::shared_ptr< Expansion2D > Expansion2DSharedPtr
Definition: Expansion1D.h:49
boost::shared_ptr< TriExp > TriExpSharedPtr
Definition: TriExp.h:283
Nektar::MultiRegions::ExpList2D::ExpList2D ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::CompositeMap domain,
const SpatialDomains::MeshGraphSharedPtr graph3D,
const std::string  variable = "DefaultVar" 
)

Specialised constructor for Neumann boundary conditions in DisContField3D and ContField3D.

Fills the list of local expansions with the segments from the 3D mesh specified by domain. This CompositeMap contains a list of Composites which define the Neumann boundary.

See also
ExpList2D::ExpList2D(SpatialDomains::MeshGraph2D&) for details.
Parameters
domainA domain, comprising of one or more composite regions.
graph3DA mesh, containing information about the domain and the spectral/hp element expansions.

Definition at line 778 of file ExpList2D.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::MultiRegions::ExpList::CreateCollections(), Nektar::MultiRegions::e2D, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eNodalTriElec, Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::LibUtilities::BasisKey::GetNumModes(), Nektar::LibUtilities::BasisKey::GetNumPoints(), Nektar::MultiRegions::ExpList::m_coeffs, Nektar::MultiRegions::ExpList::m_globalOptParam, Nektar::MultiRegions::ExpList::m_ncoeffs, Nektar::MultiRegions::ExpList::m_npoints, Nektar::MultiRegions::ExpList::m_phys, Nektar::MultiRegions::ExpList::ReadGlobalOptimizationParameters(), SetCoeffPhysOffsets(), and Nektar::MultiRegions::ExpList::SetExpType().

782  :ExpList(pSession, graph3D)
783  {
784 
785  SetExpType(e2D);
786 
787  ASSERTL0(boost::dynamic_pointer_cast<
788  SpatialDomains::MeshGraph3D>(graph3D),
789  "Expected a MeshGraph3D object.");
790 
791  int j, elmtid=0;
792  int nel = 0;
793 
796  SpatialDomains::QuadGeomSharedPtr QuadrilateralGeom;
797 
802 
803  SpatialDomains::CompositeMap::const_iterator compIt;
804  for (compIt = domain.begin(); compIt != domain.end(); ++compIt)
805  {
806  nel += (compIt->second)->size();
807  }
808 
809  for (compIt = domain.begin(); compIt != domain.end(); ++compIt)
810  {
811  for (j = 0; j < compIt->second->size(); ++j)
812  {
813  if ((TriangleGeom = boost::dynamic_pointer_cast<
814  SpatialDomains::TriGeom>((*compIt->second)[j])))
815  {
817  = boost::dynamic_pointer_cast<
818  SpatialDomains::MeshGraph3D>(graph3D)->
819  GetFaceBasisKey(TriangleGeom, 0, variable);
821  = boost::dynamic_pointer_cast<
822  SpatialDomains::MeshGraph3D>(graph3D)->
823  GetFaceBasisKey(TriangleGeom,1,variable);
824 
825  if (graph3D->GetExpansions().begin()->second->
826  m_basisKeyVector[0].GetBasisType() ==
828  {
829  ASSERTL0(false,"This method needs sorting");
831 
833  ::AllocateSharedPtr(TriBa,TriBb,TriNb,
834  TriangleGeom);
835  Ntri->SetElmtId(elmtid++);
836  (*m_exp).push_back(Ntri);
837  }
838  else
839  {
841  ::AllocateSharedPtr(TriBa, TriBb,
842  TriangleGeom);
843  tri->SetElmtId(elmtid++);
844  (*m_exp).push_back(tri);
845  }
846 
847  m_ncoeffs
848  += (TriBa.GetNumModes()*(TriBa.GetNumModes()+1))/2
849  + TriBa.GetNumModes()*(TriBb.GetNumModes()
850  -TriBa.GetNumModes());
851  m_npoints += TriBa.GetNumPoints()*TriBb.GetNumPoints();
852  }
853  else if ((QuadrilateralGeom = boost::dynamic_pointer_cast<
854  SpatialDomains::QuadGeom>((*compIt->second)[j])))
855  {
857  = boost::dynamic_pointer_cast<
858  SpatialDomains::MeshGraph3D>(graph3D)->
859  GetFaceBasisKey(QuadrilateralGeom, 0,
860  variable);
862  = boost::dynamic_pointer_cast<
863  SpatialDomains::MeshGraph3D>(graph3D)->
864  GetFaceBasisKey(QuadrilateralGeom, 1,
865  variable);
866 
868  ::AllocateSharedPtr(QuadBa, QuadBb,
869  QuadrilateralGeom);
870  quad->SetElmtId(elmtid++);
871  (*m_exp).push_back(quad);
872 
873  m_ncoeffs += QuadBa.GetNumModes()*QuadBb.GetNumModes();
874  m_npoints += QuadBa.GetNumPoints()
875  * QuadBb.GetNumPoints();
876  }
877  else
878  {
879  ASSERTL0(false,
880  "dynamic cast to a proper Geometry2D failed");
881  }
882  }
883 
884  }
885 
886  // Setup Default optimisation information.
887  nel = GetExpSize();
890 
891  // Set up m_coeffs, m_phys and offset arrays.
895 
898  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
ExpList()
The default constructor.
Definition: ExpList.cpp:93
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:971
boost::shared_ptr< QuadGeom > QuadGeomSharedPtr
Definition: QuadGeom.h:54
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:926
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:909
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1837
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:887
int GetNumPoints() const
Return points order at which basis is defined.
Definition: Basis.h:128
boost::shared_ptr< QuadExp > QuadExpSharedPtr
Definition: QuadExp.h:284
void SetCoeffPhysOffsets(void)
Definition of the total number of degrees of freedom and quadrature points and offsets to access data...
Definition: ExpList2D.cpp:1096
boost::shared_ptr< GeometryVector > Composite
Definition: MeshGraph.h:111
boost::shared_ptr< TriGeom > TriGeomSharedPtr
Definition: TriGeom.h:58
void CreateCollections(Collections::ImplementationType ImpType=Collections::eNoImpType)
Construct collections of elements containing a single element type and polynomial order from the list...
Definition: ExpList.cpp:2723
int GetNumModes() const
Returns the order of the basis.
Definition: Basis.h:84
Lagrange for SEM basis .
Definition: BasisType.h:53
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:211
Describes the specification for a Basis.
Definition: Basis.h:50
2D Nodal Electrostatic Points on a Triangle
Definition: PointsType.h:68
boost::shared_ptr< TriExp > TriExpSharedPtr
Definition: TriExp.h:283
boost::shared_ptr< NodalTriExp > NodalTriExpSharedPtr
Definition: NodalTriExp.h:379
Nektar::MultiRegions::ExpList2D::~ExpList2D ( )
virtual

Destructor.

Definition at line 78 of file ExpList2D.cpp.

79  {
80  }

Member Function Documentation

void Nektar::MultiRegions::ExpList2D::SetCoeffPhysOffsets ( void  )
private

Definition of the total number of degrees of freedom and quadrature points and offsets to access datax.

Each expansion (local element) is processed in turn to determine the number of coefficients and physical data points it contributes to the domain. Three arrays, m_coeff_offset, m_phys_offset and m_offset_elmt_id, are initialised and updated to store the data offsets of each element in the m_coeffs and m_phys arrays, and the element id that each consecutive block is associated respectively.

Definition at line 1096 of file ExpList2D.cpp.

References Nektar::LibUtilities::eQuadrilateral, Nektar::LibUtilities::eTriangle, Nektar::MultiRegions::ExpList::m_coeff_offset, Nektar::MultiRegions::ExpList::m_exp, Nektar::MultiRegions::ExpList::m_ncoeffs, Nektar::MultiRegions::ExpList::m_npoints, Nektar::MultiRegions::ExpList::m_offset_elmt_id, and Nektar::MultiRegions::ExpList::m_phys_offset.

Referenced by ExpList2D().

1097  {
1098  int i;
1099 
1100  // Set up offset information and array sizes
1102  m_phys_offset = Array<OneD,int>(m_exp->size());
1104 
1105  m_ncoeffs = m_npoints = 0;
1106 
1107  int cnt = 0;
1108  for(i = 0; i < m_exp->size(); ++i)
1109  {
1110  if((*m_exp)[i]->DetShapeType() == LibUtilities::eTriangle)
1111  {
1113  m_phys_offset [i] = m_npoints;
1114  m_offset_elmt_id[cnt++] = i;
1115  m_ncoeffs += (*m_exp)[i]->GetNcoeffs();
1116  m_npoints += (*m_exp)[i]->GetTotPoints();
1117  }
1118  }
1119 
1120  for(i = 0; i < m_exp->size(); ++i)
1121  {
1122  if((*m_exp)[i]->DetShapeType() == LibUtilities::eQuadrilateral)
1123  {
1125  m_phys_offset [i] = m_npoints;
1126  m_offset_elmt_id[cnt++] = i;
1127  m_ncoeffs += (*m_exp)[i]->GetNcoeffs();
1128  m_npoints += (*m_exp)[i]->GetTotPoints();
1129  }
1130  }
1131  }
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:958
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:947
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:961
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:887
Array< OneD, int > m_offset_elmt_id
Array containing the element id m_offset_elmt_id[n] that the n^th consecutive block of data in m_coef...
Definition: ExpList.h:969
void Nektar::MultiRegions::ExpList2D::v_GetNormals ( Array< OneD, Array< OneD, NekDouble > > &  normals)
protectedvirtual

For each local element, copy the normals stored in the element list into the array normals.

Parameters
normalsMultidimensional array in which to copy normals to. Must have dimension equal to or larger than the spatial dimension of the elements.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1014 of file ExpList2D.cpp.

References Nektar::MultiRegions::AlignFace(), ASSERTL1, Nektar::StdRegions::eDir1FwdDir2_Dir2FwdDir1, Nektar::MultiRegions::ExpList::GetCoordim(), Nektar::LocalRegions::Expansion2D::GetLeftAdjacentElementExp(), Nektar::LibUtilities::BasisKey::GetNumPoints(), Nektar::LibUtilities::BasisKey::GetPointsKey(), Nektar::LibUtilities::Interp2D(), Nektar::MultiRegions::ExpList::m_exp, and Nektar::MultiRegions::ExpList::m_phys_offset.

1016  {
1018  int i, j;
1019  const int coordim = GetCoordim(0);
1020 
1021  ASSERTL1(normals.num_elements() >= coordim,
1022  "Output vector does not have sufficient dimensions to "
1023  "match coordim");
1024 
1025  // Process each expansion.
1026  for (i = 0; i < m_exp->size(); ++i)
1027  {
1028  LocalRegions::Expansion2DSharedPtr traceExp = (*m_exp)[i]->as<
1031  traceExp->GetLeftAdjacentElementExp();
1032 
1033  // Get the number of points and normals for this expansion.
1034  int faceNum = traceExp->GetLeftAdjacentElementFace();
1035  int offset = m_phys_offset[i];
1036 
1037  const Array<OneD, const Array<OneD, NekDouble> > &locNormals
1038  = exp3D->GetFaceNormal(faceNum);
1039 
1040  // Project normals from 3D element onto the same orientation as
1041  // the trace expansion.
1042  StdRegions::Orientation orient = exp3D->GetForient(faceNum);
1043 
1044 
1045  int fromid0,fromid1;
1046 
1048  {
1049  fromid0 = 0;
1050  fromid1 = 1;
1051  }
1052  else
1053  {
1054  fromid0 = 1;
1055  fromid1 = 0;
1056  }
1057 
1058  LibUtilities::BasisKey faceBasis0
1059  = exp3D->DetFaceBasisKey(faceNum, fromid0);
1060  LibUtilities::BasisKey faceBasis1
1061  = exp3D->DetFaceBasisKey(faceNum, fromid1);
1062  LibUtilities::BasisKey traceBasis0
1063  = traceExp->GetBasis(0)->GetBasisKey();
1064  LibUtilities::BasisKey traceBasis1
1065  = traceExp->GetBasis(1)->GetBasisKey();
1066 
1067  const int faceNq0 = faceBasis0.GetNumPoints();
1068  const int faceNq1 = faceBasis1.GetNumPoints();
1069 
1070  for (j = 0; j < coordim; ++j)
1071  {
1072  Array<OneD, NekDouble> traceNormals(faceNq0 * faceNq1);
1073  AlignFace(orient, faceNq0, faceNq1,
1074  locNormals[j], traceNormals);
1076  faceBasis0.GetPointsKey(),
1077  faceBasis1.GetPointsKey(),
1078  traceNormals,
1079  traceBasis0.GetPointsKey(),
1080  traceBasis1.GetPointsKey(),
1081  tmp = normals[j]+offset);
1082  }
1083  }
1084  }
void Interp2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
this function interpolates a 2D function evaluated at the quadrature points of the 2D basis...
Definition: Interp.cpp:116
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:947
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:961
boost::shared_ptr< Expansion3D > Expansion3DSharedPtr
Definition: Expansion2D.h:48
int GetNumPoints() const
Return points order at which basis is defined.
Definition: Basis.h:128
Expansion3DSharedPtr GetLeftAdjacentElementExp() const
Definition: Expansion2D.h:193
PointsKey GetPointsKey() const
Return distribution of points.
Definition: Basis.h:145
int GetCoordim(int eid)
This function returns the dimension of the coordinates of the element eid.
Definition: ExpList.h:1735
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
Describes the specification for a Basis.
Definition: Basis.h:50
boost::shared_ptr< Expansion2D > Expansion2DSharedPtr
Definition: Expansion1D.h:49
void AlignFace(const StdRegions::Orientation orient, const int nquad1, const int nquad2, const Array< OneD, const NekDouble > &in, Array< OneD, NekDouble > &out)
Helper function to re-align face to a given orientation.
Definition: ExpList2D.cpp:942
void Nektar::MultiRegions::ExpList2D::v_PhysGalerkinProjection1DScaled ( const NekDouble  scale,
const Array< OneD, NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
privatevirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1274 of file ExpList2D.cpp.

References Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::MultiRegions::ExpList::m_exp, and Nektar::LibUtilities::PhysGalerkinProject2D().

1278  {
1279  int cnt,cnt1;
1280 
1281  cnt = cnt1 = 0;
1282  for(int i = 0; i < GetExpSize(); ++i)
1283  {
1284  // get new points key
1285  int pt0 = (*m_exp)[i]->GetNumPoints(0);
1286  int pt1 = (*m_exp)[i]->GetNumPoints(1);
1287  int npt0 = (int) pt0*scale;
1288  int npt1 = (int) pt1*scale;
1289 
1290  LibUtilities::PointsKey newPointsKey0(npt0,
1291  (*m_exp)[i]->GetPointsType(0));
1292  LibUtilities::PointsKey newPointsKey1(npt1,
1293  (*m_exp)[i]->GetPointsType(1));
1294 
1295  // Project points;
1297  newPointsKey0,
1298  newPointsKey1,
1299  &inarray[cnt],
1300  (*m_exp)[i]->GetBasis(0)->GetPointsKey(),
1301  (*m_exp)[i]->GetBasis(1)->GetPointsKey(),
1302  &outarray[cnt1]);
1303 
1304  cnt += npt0*npt1;
1305  cnt1 += pt0*pt1;
1306  }
1307 
1308  }
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1837
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:947
Defines a specification for a set of points.
Definition: Points.h:58
void PhysGalerkinProject2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
void Nektar::MultiRegions::ExpList2D::v_PhysInterp1DScaled ( const NekDouble  scale,
const Array< OneD, NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
privatevirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1242 of file ExpList2D.cpp.

References Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::LibUtilities::Interp2D(), and Nektar::MultiRegions::ExpList::m_exp.

1246  {
1247  int cnt,cnt1;
1248 
1249  cnt = cnt1 = 0;
1250  for(int i = 0; i < GetExpSize(); ++i)
1251  {
1252  // get new points key
1253  int pt0 = (*m_exp)[i]->GetNumPoints(0);
1254  int pt1 = (*m_exp)[i]->GetNumPoints(1);
1255  int npt0 = (int) pt0*scale;
1256  int npt1 = (int) pt1*scale;
1257 
1258  LibUtilities::PointsKey newPointsKey0(npt0,
1259  (*m_exp)[i]->GetPointsType(0));
1260  LibUtilities::PointsKey newPointsKey1(npt1,
1261  (*m_exp)[i]->GetPointsType(1));
1262 
1263  // Interpolate points;
1264  LibUtilities::Interp2D((*m_exp)[i]->GetBasis(0)->GetPointsKey(),
1265  (*m_exp)[i]->GetBasis(1)->GetPointsKey(),
1266  &inarray[cnt],newPointsKey0,
1267  newPointsKey1,&outarray[cnt1]);
1268 
1269  cnt += pt0*pt1;
1270  cnt1 += npt0*npt1;
1271  }
1272  }
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1837
void Interp2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
this function interpolates a 2D function evaluated at the quadrature points of the 2D basis...
Definition: Interp.cpp:116
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:947
Defines a specification for a set of points.
Definition: Points.h:58
void Nektar::MultiRegions::ExpList2D::v_ReadGlobalOptimizationParameters ( )
privatevirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1150 of file ExpList2D.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::LibUtilities::eTriangle, Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::MultiRegions::ExpList::m_exp, Nektar::MultiRegions::ExpList::m_globalOptParam, and Nektar::MultiRegions::ExpList::m_session.

1151  {
1152  Array<OneD, int> NumShape(2,0);
1153 
1154  for(int i = 0; i < GetExpSize(); ++i)
1155  {
1156  if((*m_exp)[i]->DetShapeType() == LibUtilities::eTriangle)
1157  {
1158  NumShape[0] += 1;
1159  }
1160  else // Quadrilateral element
1161  {
1162  NumShape[1] += 1;
1163  }
1164  }
1165 
1167  ::AllocateSharedPtr(m_session,2,NumShape);
1168  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:971
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1837
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:947
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:880
void Nektar::MultiRegions::ExpList2D::v_SetUpPhysNormals ( )
privatevirtual

Set up the normals on each expansion.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1137 of file ExpList2D.cpp.

References Nektar::MultiRegions::ExpList::m_exp.

1138  {
1139  int i, j;
1140  for (i = 0; i < m_exp->size(); ++i)
1141  {
1142  for (j = 0; j < (*m_exp)[i]->GetNedges(); ++j)
1143  {
1144  (*m_exp)[i]->ComputeEdgeNormal(j);
1145  }
1146  }
1147  }
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:947
void Nektar::MultiRegions::ExpList2D::v_Upwind ( const Array< OneD, const NekDouble > &  Vn,
const Array< OneD, const NekDouble > &  Fwd,
const Array< OneD, const NekDouble > &  Bwd,
Array< OneD, NekDouble > &  Upwind 
)
protectedvirtual

Upwind the Fwd and Bwd states based on the one- dimensional normal velocity field given by Vn.

One-dimensional upwind.

Parameters
VnVelocity field.
FwdLeft state.
BwdRight state.
UpwindOutput vector.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 907 of file ExpList2D.cpp.

References Nektar::MultiRegions::ExpList::m_exp, and Nektar::MultiRegions::ExpList::m_phys_offset.

912  {
913  int i,j,f_npoints,offset;
914 
915  // Process each expansion.
916  for(i = 0; i < m_exp->size(); ++i)
917  {
918  // Get the number of points and the data offset.
919  f_npoints = (*m_exp)[i]->GetNumPoints(0)*
920  (*m_exp)[i]->GetNumPoints(1);
921  offset = m_phys_offset[i];
922 
923  // Process each point in the expansion.
924  for(j = 0; j < f_npoints; ++j)
925  {
926  // Upwind based on one-dimensional velocity.
927  if(Vn[offset + j] > 0.0)
928  {
929  Upwind[offset + j] = Fwd[offset + j];
930  }
931  else
932  {
933  Upwind[offset + j] = Bwd[offset + j];
934  }
935  }
936  }
937  }
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:947
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:961
void Nektar::MultiRegions::ExpList2D::v_WriteVtkPieceHeader ( std::ostream &  outfile,
int  expansion,
int  istrip 
)
privatevirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1170 of file ExpList2D.cpp.

1174  {
1175  int i,j;
1176  int nquad0 = (*m_exp)[expansion]->GetNumPoints(0);
1177  int nquad1 = (*m_exp)[expansion]->GetNumPoints(1);
1178  int ntot = nquad0*nquad1;
1179  int ntotminus = (nquad0-1)*(nquad1-1);
1180 
1181  Array<OneD,NekDouble> coords[3];
1182  coords[0] = Array<OneD,NekDouble>(ntot,0.0);
1183  coords[1] = Array<OneD,NekDouble>(ntot,0.0);
1184  coords[2] = Array<OneD,NekDouble>(ntot,0.0);
1185  (*m_exp)[expansion]->GetCoords(coords[0],coords[1],coords[2]);
1186 
1187  outfile << " <Piece NumberOfPoints=\""
1188  << ntot << "\" NumberOfCells=\""
1189  << ntotminus << "\">" << endl;
1190  outfile << " <Points>" << endl;
1191  outfile << " <DataArray type=\"Float64\" "
1192  << "NumberOfComponents=\"3\" format=\"ascii\">" << endl;
1193  outfile << " ";
1194  for (i = 0; i < ntot; ++i)
1195  {
1196  for (j = 0; j < 3; ++j)
1197  {
1198  outfile << setprecision(8) << scientific
1199  << (float)coords[j][i] << " ";
1200  }
1201  outfile << endl;
1202  }
1203  outfile << endl;
1204  outfile << " </DataArray>" << endl;
1205  outfile << " </Points>" << endl;
1206  outfile << " <Cells>" << endl;
1207  outfile << " <DataArray type=\"Int32\" "
1208  << "Name=\"connectivity\" format=\"ascii\">" << endl;
1209  for (i = 0; i < nquad0-1; ++i)
1210  {
1211  for (j = 0; j < nquad1-1; ++j)
1212  {
1213  outfile << j*nquad0 + i << " "
1214  << j*nquad0 + i + 1 << " "
1215  << (j+1)*nquad0 + i + 1 << " "
1216  << (j+1)*nquad0 + i << endl;
1217  }
1218  }
1219  outfile << endl;
1220  outfile << " </DataArray>" << endl;
1221  outfile << " <DataArray type=\"Int32\" "
1222  << "Name=\"offsets\" format=\"ascii\">" << endl;
1223  for (i = 0; i < ntotminus; ++i)
1224  {
1225  outfile << i*4+4 << " ";
1226  }
1227  outfile << endl;
1228  outfile << " </DataArray>" << endl;
1229  outfile << " <DataArray type=\"UInt8\" "
1230  << "Name=\"types\" format=\"ascii\">" << endl;
1231  for (i = 0; i < ntotminus; ++i)
1232  {
1233  outfile << "9 ";
1234  }
1235  outfile << endl;
1236  outfile << " </DataArray>" << endl;
1237  outfile << " </Cells>" << endl;
1238  outfile << " <PointData>" << endl;
1239  }