Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 ExpList2D &In, const std::vector< unsigned int > &eIDs, const bool DeclareCoeffPhysArrays=true, const Collections::ImplementationType ImpType=Collections::eNoImpType)
 Constructor copying only elements defined in eIds. More...
 
 ExpList2D (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &graph2D, const bool DelcareCoeffPhysArrays=true, const std::string &var="DefaultVar", const Collections::ImplementationType ImpType=Collections::eNoImpType)
 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, const Collections::ImplementationType ImpType=Collections::eNoImpType)
 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, const Collections::ImplementationType ImpType=Collections::eNoImpType)
 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", const Collections::ImplementationType ImpType=Collections::eNoImpType)
 
 ExpList2D (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::CompositeMap &domain, const SpatialDomains::MeshGraphSharedPtr &graph3D, const std::string variable="DefaultVar", const Collections::ImplementationType ImpType=Collections::eNoImpType)
 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 std::vector< unsigned int > &eIDs, const bool DeclareCoeffPhysArrays=true)
 Constructor copying only elements defined in eIds. 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, const bool PhysSpaceForcing=true)
 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 DealiasedDotProd (const Array< OneD, Array< OneD, NekDouble > > &inarray1, const Array< OneD, Array< OneD, NekDouble > > &inarray2, Array< OneD, 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 NormVectorIProductWRTBase (Array< OneD, Array< OneD, NekDouble > > &V, Array< OneD, NekDouble > &outarray)
 
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 FillBndCondFromField (const int nreg)
 Fill Bnd Condition expansion in nreg from the values stored in expansion. More...
 
void LocalToGlobal (bool useComm=true)
 Gathers the global coefficients $\boldsymbol{\hat{u}}_g$ from the local coefficients $\boldsymbol{\hat{u}}_l$. More...
 
void LocalToGlobal (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool useComm=true)
 
void GlobalToLocal (void)
 Scatters from the global coefficients $\boldsymbol{\hat{u}}_g$ to the local coefficients $\boldsymbol{\hat{u}}_l$. More...
 
void GlobalToLocal (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
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::ExpansionVector
GetExp () 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...
 
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)
 
void CurlCurl (Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
 
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)
 
const std::vector< bool > & GetLeftAdjacentFaces (void) const
 
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 GetBndElmtExpansion (int i, boost::shared_ptr< ExpList > &result, const bool DeclareCoeffPhysArrays=true)
 
void ExtractElmtToBndPhys (int i, Array< OneD, NekDouble > &elmt, Array< OneD, NekDouble > &boundary)
 
void ExtractPhysToBndElmt (int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bndElmt)
 
void ExtractPhysToBnd (int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bnd)
 
void GetBoundaryNormals (int i, Array< OneD, Array< OneD, NekDouble > > &normals)
 
void GeneralGetFieldDefinitions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef, int NumHomoDir=0, Array< OneD, LibUtilities::BasisSharedPtr > &HomoBasis=LibUtilities::NullBasisSharedPtr1DArray, std::vector< NekDouble > &HomoLen=LibUtilities::NullNekDoubleVector, bool homoStrips=false, std::vector< unsigned int > &HomoSIDs=LibUtilities::NullUnsignedIntVector, std::vector< unsigned int > &HomoZIDs=LibUtilities::NullUnsignedIntVector, std::vector< unsigned int > &HomoYIDs=LibUtilities::NullUnsignedIntVector)
 
const
NekOptimize::GlobalOptParamSharedPtr
GetGlobalOptParam (void)
 
std::map< int,
RobinBCInfoSharedPtr
GetRobinBCInfo ()
 
void GetPeriodicEntities (PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces=NullPeriodicMap)
 
std::vector
< LibUtilities::FieldDefinitionsSharedPtr
GetFieldDefinitions ()
 
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::SessionReader
GetSession () const
 Returns the session object. More...
 
boost::shared_ptr
< LibUtilities::Comm
GetComm ()
 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...
 
void ClearGlobalLinSysManager (void)
 

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
void SetCoeffPhysOffsets ()
 Definition of the total number of degrees of freedom and quadrature points and offsets to access data. More...
 
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 const std::vector< bool > & v_GetLeftAdjacentFaces (void) const
 
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, const bool PhysSpaceForcing)
 
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_FillBndCondFromField (const int nreg)
 
virtual void v_Reset ()
 Reset geometry information, metrics, matrix managers and geometry information. More...
 
virtual void v_LocalToGlobal (bool UseComm)
 
virtual void v_LocalToGlobal (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool UseComm)
 
virtual void v_GlobalToLocal (void)
 
virtual void v_GlobalToLocal (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
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_CurlCurl (Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
 
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_DealiasedDotProd (const Array< OneD, Array< OneD, NekDouble > > &inarray1, const Array< OneD, Array< OneD, NekDouble > > &inarray2, Array< OneD, 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_NormVectorIProductWRTBase (Array< OneD, Array< OneD, NekDouble > > &V, Array< OneD, NekDouble > &outarray)
 
virtual void v_GetBoundaryToElmtMap (Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
 
virtual void v_GetBndElmtExpansion (int i, boost::shared_ptr< ExpList > &result, const bool DeclareCoeffPhysArrays)
 
virtual void v_ExtractElmtToBndPhys (int i, Array< OneD, NekDouble > &elmt, Array< OneD, NekDouble > &boundary)
 
virtual void v_ExtractPhysToBndElmt (int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bndElmt)
 
virtual void v_ExtractPhysToBnd (int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bnd)
 
virtual void v_GetBoundaryNormals (int i, Array< OneD, Array< OneD, NekDouble > > &normals)
 
virtual std::vector
< LibUtilities::FieldDefinitionsSharedPtr
v_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
NekDouble
v_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)
 
virtual void v_ClearGlobalLinSysManager (void)
 
void ExtractFileBCs (const std::string &fileName, LibUtilities::CommSharedPtr comm, const std::string &varName, const boost::shared_ptr< ExpList > locExpList)
 

Private Member Functions

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::ExpansionVector
m_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...
 
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
 
BlockMatrixMapShPtr m_blockMat
 
bool m_WaveSpace
 
boost::unordered_map< int, int > m_elmtToExpId
 Mapping from geometry ID of element to index inside m_exp. More...
 

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 69 of file ExpList2D.cpp.

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

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

Copy constructor.

Parameters
InExpList2D object to copy.

Definition at line 87 of file ExpList2D.cpp.

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

89  :
90  ExpList(In,DeclareCoeffPhysArrays)
91  {
92  SetExpType(e2D);
93  }
ExpList()
The default constructor.
Definition: ExpList.cpp:95
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:277
Nektar::MultiRegions::ExpList2D::ExpList2D ( const ExpList2D In,
const std::vector< unsigned int > &  eIDs,
const bool  DeclareCoeffPhysArrays = true,
const Collections::ImplementationType  ImpType = Collections::eNoImpType 
)

Constructor copying only elements defined in eIds.

Parameters
InExpList2D object to copy.
eIDsId of elements that should be copied.

Definition at line 99 of file ExpList2D.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::MultiRegions::ExpList::CreateCollections(), Nektar::MultiRegions::e2D, Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::MultiRegions::ExpList::m_globalOptParam, Nektar::MultiRegions::ExpList::ReadGlobalOptimizationParameters(), Nektar::MultiRegions::ExpList::SetCoeffPhysOffsets(), and Nektar::MultiRegions::ExpList::SetExpType().

103  :
104  ExpList(In,eIDs,DeclareCoeffPhysArrays)
105  {
106  SetExpType(e2D);
107 
108  // Setup Default optimisation information.
109  int nel = GetExpSize();
112 
113  // set up offset arrays.
115 
117  CreateCollections(ImpType);
118  }
ExpList()
The default constructor.
Definition: ExpList.cpp:95
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:1052
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2046
void SetCoeffPhysOffsets()
Definition of the total number of degrees of freedom and quadrature points and offsets to access data...
Definition: ExpList.cpp:247
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:3151
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:277
Nektar::MultiRegions::ExpList2D::ExpList2D ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr graph2D,
const bool  DeclareCoeffPhysArrays = true,
const std::string &  var = "DefaultVar",
const Collections::ImplementationType  ImpType = Collections::eNoImpType 
)

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 133 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(), Nektar::MultiRegions::ExpList::SetCoeffPhysOffsets(), and Nektar::MultiRegions::ExpList::SetExpType().

138  :
139  ExpList(pSession,graph2D)
140  {
141  SetExpType(e2D);
142 
143  int elmtid=0;
149 
150  const SpatialDomains::ExpansionMap &expansions
151  = graph2D->GetExpansions(var);
152 
153  SpatialDomains::ExpansionMap::const_iterator expIt;
154  for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
155  {
157  SpatialDomains::QuadGeomSharedPtr QuadrilateralGeom;
158 
159  if ((TriangleGeom = boost::dynamic_pointer_cast<SpatialDomains
160  ::TriGeom>(expIt->second->m_geomShPtr)))
161  {
162  LibUtilities::BasisKey TriBa
163  = expIt->second->m_basisKeyVector[0];
164  LibUtilities::BasisKey TriBb
165  = expIt->second->m_basisKeyVector[1];
166 
167  // This is not elegantly implemented needs re-thinking.
168  if (TriBa.GetBasisType() == LibUtilities::eGLL_Lagrange)
169  {
170  LibUtilities::BasisKey newBa(LibUtilities::eOrtho_A,
171  TriBa.GetNumModes(),
172  TriBa.GetPointsKey());
173 
176  ::AllocateSharedPtr(newBa,TriBb,TriNb,
177  TriangleGeom);
178  Ntri->SetElmtId(elmtid++);
179  (*m_exp).push_back(Ntri);
180  }
181  else
182  {
184  ::AllocateSharedPtr(TriBa,TriBb,
185  TriangleGeom);
186  tri->SetElmtId(elmtid++);
187  (*m_exp).push_back(tri);
188  }
189  m_ncoeffs += (TriBa.GetNumModes()*(TriBa.GetNumModes()+1))/2
190  + TriBa.GetNumModes()*(TriBb.GetNumModes()
191  -TriBa.GetNumModes());
192  m_npoints += TriBa.GetNumPoints()*TriBb.GetNumPoints();
193  }
194  else if ((QuadrilateralGeom = boost::dynamic_pointer_cast<
195  SpatialDomains::QuadGeom>(expIt->second->m_geomShPtr)))
196  {
197  LibUtilities::BasisKey QuadBa
198  = expIt->second->m_basisKeyVector[0];
199  LibUtilities::BasisKey QuadBb
200  = expIt->second->m_basisKeyVector[1];
201 
203  ::AllocateSharedPtr(QuadBa,QuadBb,
204  QuadrilateralGeom);
205  quad->SetElmtId(elmtid++);
206  (*m_exp).push_back(quad);
207 
208  m_ncoeffs += QuadBa.GetNumModes()*QuadBb.GetNumModes();
209  m_npoints += QuadBa.GetNumPoints()*QuadBb.GetNumPoints();
210  }
211  else
212  {
213  ASSERTL0(false, "dynamic cast to a proper Geometry2D "
214  "failed");
215  }
216 
217  }
218 
219  // set up element numbering
220  for(int i = 0; i < (*m_exp).size(); ++i)
221  {
222  (*m_exp)[i]->SetElmtId(i);
223  }
224 
225  // Setup Default optimisation information.
226  int nel = GetExpSize();
229 
230 
231  // set up offset arrays.
233 
234  if (DeclareCoeffPhysArrays)
235  {
236  // Set up m_coeffs, m_phys.
237  m_coeffs = Array<OneD, NekDouble>(m_ncoeffs);
238  m_phys = Array<OneD, NekDouble>(m_npoints);
239  }
240 
242  CreateCollections(ImpType);
243  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
ExpList()
The default constructor.
Definition: ExpList.cpp:95
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:1052
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:1015
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:998
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2046
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:976
void SetCoeffPhysOffsets()
Definition of the total number of degrees of freedom and quadrature points and offsets to access data...
Definition: ExpList.cpp:247
Principle Orthogonal Functions .
Definition: BasisType.h:46
boost::shared_ptr< QuadExp > QuadExpSharedPtr
Definition: QuadExp.h:297
boost::shared_ptr< GeometryVector > Composite
Definition: MeshGraph.h:114
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:3151
Lagrange for SEM basis .
Definition: BasisType.h:53
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:277
2D Nodal Electrostatic Points on a Triangle
Definition: PointsType.h:70
std::map< int, ExpansionShPtr > ExpansionMap
Definition: MeshGraph.h:174
boost::shared_ptr< TriExp > TriExpSharedPtr
Definition: TriExp.h:295
boost::shared_ptr< NodalTriExp > NodalTriExpSharedPtr
Definition: NodalTriExp.h:381
Nektar::MultiRegions::ExpList2D::ExpList2D ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::ExpansionMap expansions,
const bool  DeclareCoeffPhysArrays = true,
const Collections::ImplementationType  ImpType = Collections::eNoImpType 
)

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 260 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(), Nektar::MultiRegions::ExpList::SetCoeffPhysOffsets(), and Nektar::MultiRegions::ExpList::SetExpType().

264  :
265  ExpList(pSession)
266  {
267  SetExpType(e2D);
268 
269  int elmtid=0;
275 
277  for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
278  {
280  SpatialDomains::QuadGeomSharedPtr QuadrilateralGeom;
281 
282  if ((TriangleGeom = boost::dynamic_pointer_cast<SpatialDomains
283  ::TriGeom>(expIt->second->m_geomShPtr)))
284  {
285  LibUtilities::BasisKey TriBa
286  = expIt->second->m_basisKeyVector[0];
287  LibUtilities::BasisKey TriBb
288  = expIt->second->m_basisKeyVector[1];
289 
290  // This is not elegantly implemented needs re-thinking.
291  if (TriBa.GetBasisType() == LibUtilities::eGLL_Lagrange)
292  {
293  LibUtilities::BasisKey newBa(LibUtilities::eOrtho_A,
294  TriBa.GetNumModes(),
295  TriBa.GetPointsKey());
296 
299  ::AllocateSharedPtr(newBa,TriBb,TriNb,
300  TriangleGeom);
301  Ntri->SetElmtId(elmtid++);
302  (*m_exp).push_back(Ntri);
303  }
304  else
305  {
307  ::AllocateSharedPtr(TriBa,TriBb,
308  TriangleGeom);
309  tri->SetElmtId(elmtid++);
310  (*m_exp).push_back(tri);
311  }
312  m_ncoeffs += (TriBa.GetNumModes()*(TriBa.GetNumModes()+1))/2
313  + TriBa.GetNumModes()*(TriBb.GetNumModes()
314  -TriBa.GetNumModes());
315  m_npoints += TriBa.GetNumPoints()*TriBb.GetNumPoints();
316  }
317  else if ((QuadrilateralGeom = boost::dynamic_pointer_cast<
318  SpatialDomains::QuadGeom>(expIt->second->m_geomShPtr)))
319  {
320  LibUtilities::BasisKey QuadBa
321  = expIt->second->m_basisKeyVector[0];
322  LibUtilities::BasisKey QuadBb
323  = expIt->second->m_basisKeyVector[1];
324 
326  ::AllocateSharedPtr(QuadBa,QuadBb,
327  QuadrilateralGeom);
328  quad->SetElmtId(elmtid++);
329  (*m_exp).push_back(quad);
330 
331  m_ncoeffs += QuadBa.GetNumModes()*QuadBb.GetNumModes();
332  m_npoints += QuadBa.GetNumPoints()*QuadBb.GetNumPoints();
333  }
334  else
335  {
336  ASSERTL0(false, "dynamic cast to a proper Geometry2D "
337  "failed");
338  }
339 
340  }
341 
342  // Setup Default optimisation information.
343  int nel = GetExpSize();
346 
347 
348  // set up offset arrays.
350 
351  if (DeclareCoeffPhysArrays)
352  {
353  // Set up m_coeffs, m_phys.
354  m_coeffs = Array<OneD, NekDouble>(m_ncoeffs);
355  m_phys = Array<OneD, NekDouble>(m_npoints);
356  }
357 
359  CreateCollections(ImpType);
360  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
std::map< int, ExpansionShPtr >::const_iterator ExpansionMapConstIter
Definition: MeshGraph.h:176
ExpList()
The default constructor.
Definition: ExpList.cpp:95
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:1052
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:1015
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:998
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2046
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:976
void SetCoeffPhysOffsets()
Definition of the total number of degrees of freedom and quadrature points and offsets to access data...
Definition: ExpList.cpp:247
Principle Orthogonal Functions .
Definition: BasisType.h:46
boost::shared_ptr< QuadExp > QuadExpSharedPtr
Definition: QuadExp.h:297
boost::shared_ptr< GeometryVector > Composite
Definition: MeshGraph.h:114
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:3151
Lagrange for SEM basis .
Definition: BasisType.h:53
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:277
2D Nodal Electrostatic Points on a Triangle
Definition: PointsType.h:70
boost::shared_ptr< TriExp > TriExpSharedPtr
Definition: TriExp.h:295
boost::shared_ptr< NodalTriExp > NodalTriExpSharedPtr
Definition: NodalTriExp.h:381
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,
const Collections::ImplementationType  ImpType = Collections::eNoImpType 
)

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 388 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(), Nektar::MultiRegions::ExpList::SetCoeffPhysOffsets(), Nektar::MultiRegions::ExpList::SetExpType(), and Nektar::LibUtilities::SIZE_PointsType.

396  :
397  ExpList(pSession, graph2D)
398  {
399  SetExpType(e2D);
400 
401  int elmtid=0;
406 
407  const SpatialDomains::ExpansionMap &expansions =
408  graph2D->GetExpansions();
409  m_ncoeffs = 0;
410  m_npoints = 0;
411 
412  m_physState = false;
413 
414  SpatialDomains::ExpansionMap::const_iterator expIt;
415  for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
416  {
418  SpatialDomains::QuadGeomSharedPtr QuadrilateralGeom;
419 
420  if ((TriangleGeom = boost::dynamic_pointer_cast<SpatialDomains::
421  TriGeom>(expIt->second->m_geomShPtr)))
422  {
423  if (TriNb < LibUtilities::SIZE_PointsType)
424  {
426  AllocateSharedPtr(TriBa, TriBb, TriNb,
427  TriangleGeom);
428  Ntri->SetElmtId(elmtid++);
429  (*m_exp).push_back(Ntri);
430  }
431  else
432  {
434  AllocateSharedPtr(TriBa, TriBb, TriangleGeom);
435  tri->SetElmtId(elmtid++);
436  (*m_exp).push_back(tri);
437  }
438 
439  m_ncoeffs += (TriBa.GetNumModes()*(TriBa.GetNumModes()+1))/2
440  + TriBa.GetNumModes() * (TriBb.GetNumModes() -
441  TriBa.GetNumModes());
442  m_npoints += TriBa.GetNumPoints()*TriBb.GetNumPoints();
443  }
444  else if ((QuadrilateralGeom = boost::dynamic_pointer_cast<
445  SpatialDomains::QuadGeom>(expIt->second->m_geomShPtr)))
446  {
448  AllocateSharedPtr(QuadBa, QuadBb, QuadrilateralGeom);
449  quad->SetElmtId(elmtid++);
450  (*m_exp).push_back(quad);
451 
452  m_ncoeffs += QuadBa.GetNumModes()*QuadBb.GetNumModes();
453  m_npoints += QuadBa.GetNumPoints()*QuadBb.GetNumPoints();
454  }
455  else
456  {
457  ASSERTL0(false,
458  "dynamic cast to a proper Geometry2D failed");
459  }
460 
461  }
462 
463  // Setup Default optimisation information.
464  int nel = GetExpSize();
467 
468  // Set up m_coeffs, m_phys and offset arrays.
470  m_coeffs = Array<OneD, NekDouble>(m_ncoeffs);
471  m_phys = Array<OneD, NekDouble>(m_npoints);
472 
474  CreateCollections(ImpType);
475  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
ExpList()
The default constructor.
Definition: ExpList.cpp:95
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:1052
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:1015
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:998
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2046
bool m_physState
The state of the array m_phys.
Definition: ExpList.h:1024
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:976
void SetCoeffPhysOffsets()
Definition of the total number of degrees of freedom and quadrature points and offsets to access data...
Definition: ExpList.cpp:247
boost::shared_ptr< QuadExp > QuadExpSharedPtr
Definition: QuadExp.h:297
boost::shared_ptr< GeometryVector > Composite
Definition: MeshGraph.h:114
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:3151
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:277
std::map< int, ExpansionShPtr > ExpansionMap
Definition: MeshGraph.h:174
boost::shared_ptr< TriExp > TriExpSharedPtr
Definition: TriExp.h:295
boost::shared_ptr< NodalTriExp > NodalTriExpSharedPtr
Definition: NodalTriExp.h:381
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",
const Collections::ImplementationType  ImpType = Collections::eNoImpType 
)

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 493 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, Nektar::MultiRegions::ExpList::SetCoeffPhysOffsets(), Nektar::MultiRegions::ExpList::SetExpType(), and Vmath::Vsum().

502  :
503  ExpList(pSession, graph3D)
504  {
505  SetExpType(e2D);
506 
507  int i, j, id, elmtid=0;
508  set<int> facesDone;
509 
513  LocalRegions::QuadExpSharedPtr FaceQuadExp;
517 
518  // First loop over boundary conditions to renumber
519  // Dirichlet boundaries
520  for (i = 0; i < bndCond.num_elements(); ++i)
521  {
522  if (bndCond[i]->GetBoundaryConditionType()
524  {
525  for (j = 0; j < bndConstraint[i]->GetExpSize(); ++j)
526  {
527  LibUtilities::BasisKey bkey0 = bndConstraint[i]
528  ->GetExp(j)->GetBasis(0)->GetBasisKey();
529  LibUtilities::BasisKey bkey1 = bndConstraint[i]
530  ->GetExp(j)->GetBasis(1)->GetBasisKey();
531  exp2D = bndConstraint[i]->GetExp(j)
532  ->as<LocalRegions::Expansion2D>();
533  FaceGeom = exp2D->GetGeom2D();
534 
535  //if face is a quad
536  if((FaceQuadGeom = boost::dynamic_pointer_cast<
537  SpatialDomains::QuadGeom>(FaceGeom)))
538  {
540  ::AllocateSharedPtr(bkey0, bkey1, FaceQuadGeom);
541  facesDone.insert(FaceQuadGeom->GetFid());
542  FaceQuadExp->SetElmtId(elmtid++);
543  (*m_exp).push_back(FaceQuadExp);
544  }
545  //if face is a triangle
546  else if((FaceTriGeom = boost::dynamic_pointer_cast<
547  SpatialDomains::TriGeom>(FaceGeom)))
548  {
550  ::AllocateSharedPtr(bkey0, bkey1, FaceTriGeom);
551  facesDone.insert(FaceTriGeom->GetFid());
552  FaceTriExp->SetElmtId(elmtid++);
553  (*m_exp).push_back(FaceTriExp);
554  }
555  else
556  {
557  ASSERTL0(false,"dynamic cast to a proper face geometry failed");
558  }
559  }
560  }
561  }
562 
564  pair<LibUtilities::BasisKey,
565  LibUtilities::BasisKey> > > faceOrders;
567  pair<LibUtilities::BasisKey,
568  LibUtilities::BasisKey> > >::iterator it;
569 
570  for(i = 0; i < locexp.size(); ++i)
571  {
572  exp3D = locexp[i]->as<LocalRegions::Expansion3D>();
573  for (j = 0; j < exp3D->GetNfaces(); ++j)
574  {
575  FaceGeom = exp3D->GetGeom3D()->GetFace(j);
576  id = FaceGeom->GetFid();
577 
578  if(facesDone.count(id) != 0)
579  {
580  continue;
581  }
582  it = faceOrders.find(id);
583 
584  if (it == faceOrders.end())
585  {
586  LibUtilities::BasisKey face_dir0
587  = locexp[i]->DetFaceBasisKey(j,0);
588  LibUtilities::BasisKey face_dir1
589  = locexp[i]->DetFaceBasisKey(j,1);
590 
591  faceOrders.insert(
592  std::make_pair(
593  id, std::make_pair(
594  FaceGeom,
595  std::make_pair(face_dir0, face_dir1))));
596  }
597  else // variable modes/points
598  {
599  LibUtilities::BasisKey face0 =
600  locexp[i]->DetFaceBasisKey(j,0);
601  LibUtilities::BasisKey face1 =
602  locexp[i]->DetFaceBasisKey(j,1);
603  LibUtilities::BasisKey existing0 =
604  it->second.second.first;
605  LibUtilities::BasisKey existing1 =
606  it->second.second.second;
607 
608  int np11 = face0 .GetNumPoints();
609  int np12 = face1 .GetNumPoints();
610  int np21 = existing0.GetNumPoints();
611  int np22 = existing1.GetNumPoints();
612  int nm11 = face0 .GetNumModes ();
613  int nm12 = face1 .GetNumModes ();
614  int nm21 = existing0.GetNumModes ();
615  int nm22 = existing1.GetNumModes ();
616 
617  if ((np22 >= np12 || np21 >= np11) &&
618  (nm22 >= nm12 || nm21 >= nm11))
619  {
620  continue;
621  }
622  else if((np22 < np12 || np21 < np11) &&
623  (nm22 < nm12 || nm21 < nm11))
624  {
625  it->second.second.first = face0;
626  it->second.second.second = face1;
627  }
628  else
629  {
630  ASSERTL0(false,
631  "inappropriate number of points/modes (max "
632  "num of points is not set with max order)");
633  }
634  }
635  }
636  }
637 
638  LibUtilities::CommSharedPtr vComm = pSession->GetComm();
639  int nproc = vComm->GetSize(); // number of processors
640  int facepr = vComm->GetRank(); // ID processor
641 
642  if (nproc > 1)
643  {
644  int fCnt = 0;
645 
646  // Count the number of faces on each partition
647  for(i = 0; i < locexp.size(); ++i)
648  {
649  fCnt += locexp[i]->GetNfaces();
650  }
651 
652  // Set up the offset and the array that will contain the list of
653  // face IDs, then reduce this across processors.
654  Array<OneD, int> faceCnt(nproc,0);
655  faceCnt[facepr] = fCnt;
656  vComm->AllReduce(faceCnt, LibUtilities::ReduceSum);
657 
658  int totFaceCnt = Vmath::Vsum(nproc, faceCnt, 1);
659  Array<OneD, int> fTotOffsets(nproc,0);
660 
661  for (i = 1; i < nproc; ++i)
662  {
663  fTotOffsets[i] = fTotOffsets[i-1] + faceCnt[i-1];
664  }
665 
666  // Local list of the edges per element
667 
668  Array<OneD, int> FacesTotID (totFaceCnt, 0);
669  Array<OneD, int> FacesTotNm0 (totFaceCnt, 0);
670  Array<OneD, int> FacesTotNm1 (totFaceCnt, 0);
671  Array<OneD, int> FacesTotPnts0(totFaceCnt, 0);
672  Array<OneD, int> FacesTotPnts1(totFaceCnt, 0);
673 
674  int cntr = fTotOffsets[facepr];
675 
676  for(i = 0; i < locexp.size(); ++i)
677  {
678  exp3D = locexp[i]->as<LocalRegions::Expansion3D>();
679 
680  int nfaces = locexp[i]->GetNfaces();
681 
682  for(j = 0; j < nfaces; ++j, ++cntr)
683  {
684  LibUtilities::BasisKey face_dir0
685  = locexp[i]->DetFaceBasisKey(j,0);
686  LibUtilities::BasisKey face_dir1
687  = locexp[i]->DetFaceBasisKey(j,1);
688 
689  FacesTotID[cntr] = exp3D->GetGeom3D()->GetFid(j);
690  FacesTotNm0[cntr] = face_dir0.GetNumModes ();
691  FacesTotNm1[cntr] = face_dir1.GetNumModes ();
692  FacesTotPnts0[cntr] = face_dir0.GetNumPoints();
693  FacesTotPnts1[cntr] = face_dir1.GetNumPoints();
694  }
695  }
696 
697  vComm->AllReduce(FacesTotID, LibUtilities::ReduceSum);
698  vComm->AllReduce(FacesTotNm0, LibUtilities::ReduceSum);
699  vComm->AllReduce(FacesTotNm1, LibUtilities::ReduceSum);
700  vComm->AllReduce(FacesTotPnts0, LibUtilities::ReduceSum);
701  vComm->AllReduce(FacesTotPnts1, LibUtilities::ReduceSum);
702 
703  for (i = 0; i < totFaceCnt; ++i)
704  {
705  it = faceOrders.find(FacesTotID[i]);
706 
707  if (it == faceOrders.end())
708  {
709  continue;
710  }
711 
712  LibUtilities::BasisKey existing0 =
713  it->second.second.first;
714  LibUtilities::BasisKey existing1 =
715  it->second.second.second;
716  LibUtilities::BasisKey face0(
717  existing0.GetBasisType(), FacesTotNm0[i],
718  LibUtilities::PointsKey(FacesTotPnts0[i],
719  existing0.GetPointsType()));
720  LibUtilities::BasisKey face1(
721  existing1.GetBasisType(), FacesTotNm1[i],
722  LibUtilities::PointsKey(FacesTotPnts1[i],
723  existing1.GetPointsType()));
724 
725  int np11 = face0 .GetNumPoints();
726  int np12 = face1 .GetNumPoints();
727  int np21 = existing0.GetNumPoints();
728  int np22 = existing1.GetNumPoints();
729  int nm11 = face0 .GetNumModes ();
730  int nm12 = face1 .GetNumModes ();
731  int nm21 = existing0.GetNumModes ();
732  int nm22 = existing1.GetNumModes ();
733 
734  if ((np22 >= np12 || np21 >= np11) &&
735  (nm22 >= nm12 || nm21 >= nm11))
736  {
737  continue;
738  }
739  else if((np22 < np12 || np21 < np11) &&
740  (nm22 < nm12 || nm21 < nm11))
741  {
742  it->second.second.first = face0;
743  it->second.second.second = face1;
744  }
745  else
746  {
747  ASSERTL0(false,
748  "inappropriate number of points/modes (max "
749  "num of points is not set with max order)");
750  }
751  }
752  }
753 
754  for (it = faceOrders.begin(); it != faceOrders.end(); ++it)
755  {
756  FaceGeom = it->second.first;
757 
758  if ((FaceQuadGeom = boost::dynamic_pointer_cast<
759  SpatialDomains::QuadGeom>(FaceGeom)))
760  {
762  ::AllocateSharedPtr(it->second.second.first,
763  it->second.second.second,
764  FaceQuadGeom);
765  FaceQuadExp->SetElmtId(elmtid++);
766  (*m_exp).push_back(FaceQuadExp);
767  }
768  else if ((FaceTriGeom = boost::dynamic_pointer_cast<
769  SpatialDomains::TriGeom>(FaceGeom)))
770  {
772  ::AllocateSharedPtr(it->second.second.first,
773  it->second.second.second,
774  FaceTriGeom);
775  FaceTriExp->SetElmtId(elmtid++);
776  (*m_exp).push_back(FaceTriExp);
777  }
778  }
779 
780  // Setup Default optimisation information.
781  int nel = GetExpSize();
782 
785 
786  // Set up offset information and array sizes
788 
789  // Set up m_coeffs, m_phys.
790  if(DeclareCoeffPhysArrays)
791  {
792  m_coeffs = Array<OneD, NekDouble>(m_ncoeffs);
793  m_phys = Array<OneD, NekDouble>(m_npoints);
794  }
795 
796  CreateCollections(ImpType);
797  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
ExpList()
The default constructor.
Definition: ExpList.cpp:95
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:1052
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:1015
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:998
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2046
boost::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
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:976
void SetCoeffPhysOffsets()
Definition of the total number of degrees of freedom and quadrature points and offsets to access data...
Definition: ExpList.cpp:247
boost::shared_ptr< QuadExp > QuadExpSharedPtr
Definition: QuadExp.h:297
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:3151
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:277
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:737
boost::shared_ptr< Expansion2D > Expansion2DSharedPtr
Definition: Expansion1D.h:49
boost::shared_ptr< TriExp > TriExpSharedPtr
Definition: TriExp.h:295
Nektar::MultiRegions::ExpList2D::ExpList2D ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::CompositeMap domain,
const SpatialDomains::MeshGraphSharedPtr graph3D,
const std::string  variable = "DefaultVar",
const Collections::ImplementationType  ImpType = Collections::eNoImpType 
)

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 810 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(), Nektar::MultiRegions::ExpList::SetCoeffPhysOffsets(), and Nektar::MultiRegions::ExpList::SetExpType().

815  :
816  ExpList(pSession, graph3D)
817  {
818 
819  SetExpType(e2D);
820 
821  ASSERTL0(boost::dynamic_pointer_cast<
822  SpatialDomains::MeshGraph3D>(graph3D),
823  "Expected a MeshGraph3D object.");
824 
825  int j, elmtid=0;
826  int nel = 0;
827 
830  SpatialDomains::QuadGeomSharedPtr QuadrilateralGeom;
831 
836 
837  SpatialDomains::CompositeMap::const_iterator compIt;
838  for (compIt = domain.begin(); compIt != domain.end(); ++compIt)
839  {
840  nel += (compIt->second)->size();
841  }
842 
843  for (compIt = domain.begin(); compIt != domain.end(); ++compIt)
844  {
845  for (j = 0; j < compIt->second->size(); ++j)
846  {
847  if ((TriangleGeom = boost::dynamic_pointer_cast<
848  SpatialDomains::TriGeom>((*compIt->second)[j])))
849  {
850  LibUtilities::BasisKey TriBa
851  = boost::dynamic_pointer_cast<
852  SpatialDomains::MeshGraph3D>(graph3D)->
853  GetFaceBasisKey(TriangleGeom, 0, variable);
854  LibUtilities::BasisKey TriBb
855  = boost::dynamic_pointer_cast<
856  SpatialDomains::MeshGraph3D>(graph3D)->
857  GetFaceBasisKey(TriangleGeom,1,variable);
858 
859  if (graph3D->GetExpansions().begin()->second->
860  m_basisKeyVector[0].GetBasisType() ==
862  {
863  ASSERTL0(false,"This method needs sorting");
865 
867  ::AllocateSharedPtr(TriBa,TriBb,TriNb,
868  TriangleGeom);
869  Ntri->SetElmtId(elmtid++);
870  (*m_exp).push_back(Ntri);
871  }
872  else
873  {
875  ::AllocateSharedPtr(TriBa, TriBb,
876  TriangleGeom);
877  tri->SetElmtId(elmtid++);
878  (*m_exp).push_back(tri);
879  }
880 
881  m_ncoeffs
882  += (TriBa.GetNumModes()*(TriBa.GetNumModes()+1))/2
883  + TriBa.GetNumModes()*(TriBb.GetNumModes()
884  -TriBa.GetNumModes());
885  m_npoints += TriBa.GetNumPoints()*TriBb.GetNumPoints();
886  }
887  else if ((QuadrilateralGeom = boost::dynamic_pointer_cast<
888  SpatialDomains::QuadGeom>((*compIt->second)[j])))
889  {
890  LibUtilities::BasisKey QuadBa
891  = boost::dynamic_pointer_cast<
892  SpatialDomains::MeshGraph3D>(graph3D)->
893  GetFaceBasisKey(QuadrilateralGeom, 0,
894  variable);
895  LibUtilities::BasisKey QuadBb
896  = boost::dynamic_pointer_cast<
897  SpatialDomains::MeshGraph3D>(graph3D)->
898  GetFaceBasisKey(QuadrilateralGeom, 1,
899  variable);
900 
902  ::AllocateSharedPtr(QuadBa, QuadBb,
903  QuadrilateralGeom);
904  quad->SetElmtId(elmtid++);
905  (*m_exp).push_back(quad);
906 
907  m_ncoeffs += QuadBa.GetNumModes()*QuadBb.GetNumModes();
908  m_npoints += QuadBa.GetNumPoints()
909  * QuadBb.GetNumPoints();
910  }
911  else
912  {
913  ASSERTL0(false,
914  "dynamic cast to a proper Geometry2D failed");
915  }
916  }
917 
918  }
919 
920  // Setup Default optimisation information.
921  nel = GetExpSize();
924 
925  // Set up m_coeffs, m_phys and offset arrays.
927  m_coeffs = Array<OneD, NekDouble>(m_ncoeffs);
928  m_phys = Array<OneD, NekDouble>(m_npoints);
929 
931  CreateCollections(ImpType);
932  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
ExpList()
The default constructor.
Definition: ExpList.cpp:95
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:1052
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:1015
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:998
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2046
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:976
void SetCoeffPhysOffsets()
Definition of the total number of degrees of freedom and quadrature points and offsets to access data...
Definition: ExpList.cpp:247
boost::shared_ptr< QuadExp > QuadExpSharedPtr
Definition: QuadExp.h:297
boost::shared_ptr< GeometryVector > Composite
Definition: MeshGraph.h:114
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:3151
Lagrange for SEM basis .
Definition: BasisType.h:53
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:277
2D Nodal Electrostatic Points on a Triangle
Definition: PointsType.h:70
boost::shared_ptr< TriExp > TriExpSharedPtr
Definition: TriExp.h:295
boost::shared_ptr< NodalTriExp > NodalTriExpSharedPtr
Definition: NodalTriExp.h:381
Nektar::MultiRegions::ExpList2D::~ExpList2D ( )
virtual

Destructor.

Definition at line 79 of file ExpList2D.cpp.

80  {
81  }

Member Function Documentation

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 1048 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.

1050  {
1051  Array<OneD, NekDouble> tmp;
1052  int i, j;
1053  const int coordim = GetCoordim(0);
1054 
1055  ASSERTL1(normals.num_elements() >= coordim,
1056  "Output vector does not have sufficient dimensions to "
1057  "match coordim");
1058 
1059  // Process each expansion.
1060  for (i = 0; i < m_exp->size(); ++i)
1061  {
1062  LocalRegions::Expansion2DSharedPtr traceExp = (*m_exp)[i]->as<
1063  LocalRegions::Expansion2D>();
1065  traceExp->GetLeftAdjacentElementExp();
1066 
1067  // Get the number of points and normals for this expansion.
1068  int faceNum = traceExp->GetLeftAdjacentElementFace();
1069  int offset = m_phys_offset[i];
1070 
1071  const Array<OneD, const Array<OneD, NekDouble> > &locNormals
1072  = exp3D->GetFaceNormal(faceNum);
1073 
1074  // Project normals from 3D element onto the same orientation as
1075  // the trace expansion.
1076  StdRegions::Orientation orient = exp3D->GetForient(faceNum);
1077 
1078 
1079  int fromid0,fromid1;
1080 
1082  {
1083  fromid0 = 0;
1084  fromid1 = 1;
1085  }
1086  else
1087  {
1088  fromid0 = 1;
1089  fromid1 = 0;
1090  }
1091 
1092  LibUtilities::BasisKey faceBasis0
1093  = exp3D->DetFaceBasisKey(faceNum, fromid0);
1094  LibUtilities::BasisKey faceBasis1
1095  = exp3D->DetFaceBasisKey(faceNum, fromid1);
1096  LibUtilities::BasisKey traceBasis0
1097  = traceExp->GetBasis(0)->GetBasisKey();
1098  LibUtilities::BasisKey traceBasis1
1099  = traceExp->GetBasis(1)->GetBasisKey();
1100 
1101  const int faceNq0 = faceBasis0.GetNumPoints();
1102  const int faceNq1 = faceBasis1.GetNumPoints();
1103 
1104  for (j = 0; j < coordim; ++j)
1105  {
1106  Array<OneD, NekDouble> traceNormals(faceNq0 * faceNq1);
1107  AlignFace(orient, faceNq0, faceNq1,
1108  locNormals[j], traceNormals);
1110  faceBasis0.GetPointsKey(),
1111  faceBasis1.GetPointsKey(),
1112  traceNormals,
1113  traceBasis0.GetPointsKey(),
1114  traceBasis1.GetPointsKey(),
1115  tmp = normals[j]+offset);
1116  }
1117  }
1118  }
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:1036
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:1050
boost::shared_ptr< Expansion3D > Expansion3DSharedPtr
Definition: Expansion2D.h:48
int GetCoordim(int eid)
This function returns the dimension of the coordinates of the element eid.
Definition: ExpList.h:1898
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
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:976
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 1260 of file ExpList2D.cpp.

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

1264  {
1265  int cnt,cnt1;
1266 
1267  cnt = cnt1 = 0;
1268  for(int i = 0; i < GetExpSize(); ++i)
1269  {
1270  // get new points key
1271  int pt0 = (*m_exp)[i]->GetNumPoints(0);
1272  int pt1 = (*m_exp)[i]->GetNumPoints(1);
1273  int npt0 = (int) pt0*scale;
1274  int npt1 = (int) pt1*scale;
1275 
1276  LibUtilities::PointsKey newPointsKey0(npt0,
1277  (*m_exp)[i]->GetPointsType(0));
1278  LibUtilities::PointsKey newPointsKey1(npt1,
1279  (*m_exp)[i]->GetPointsType(1));
1280 
1281  // Project points;
1283  newPointsKey0,
1284  newPointsKey1,
1285  &inarray[cnt],
1286  (*m_exp)[i]->GetBasis(0)->GetPointsKey(),
1287  (*m_exp)[i]->GetBasis(1)->GetPointsKey(),
1288  &outarray[cnt1]);
1289 
1290  cnt += npt0*npt1;
1291  cnt1 += pt0*pt1;
1292  }
1293 
1294  }
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2046
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1036
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 1228 of file ExpList2D.cpp.

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

1232  {
1233  int cnt,cnt1;
1234 
1235  cnt = cnt1 = 0;
1236  for(int i = 0; i < GetExpSize(); ++i)
1237  {
1238  // get new points key
1239  int pt0 = (*m_exp)[i]->GetNumPoints(0);
1240  int pt1 = (*m_exp)[i]->GetNumPoints(1);
1241  int npt0 = (int) pt0*scale;
1242  int npt1 = (int) pt1*scale;
1243 
1244  LibUtilities::PointsKey newPointsKey0(npt0,
1245  (*m_exp)[i]->GetPointsType(0));
1246  LibUtilities::PointsKey newPointsKey1(npt1,
1247  (*m_exp)[i]->GetPointsType(1));
1248 
1249  // Interpolate points;
1250  LibUtilities::Interp2D((*m_exp)[i]->GetBasis(0)->GetPointsKey(),
1251  (*m_exp)[i]->GetBasis(1)->GetPointsKey(),
1252  &inarray[cnt],newPointsKey0,
1253  newPointsKey1,&outarray[cnt1]);
1254 
1255  cnt += pt0*pt1;
1256  cnt1 += npt0*npt1;
1257  }
1258  }
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2046
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:1036
void Nektar::MultiRegions::ExpList2D::v_ReadGlobalOptimizationParameters ( )
privatevirtual

Reimplemented from Nektar::MultiRegions::ExpList.

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

1137  {
1138  Array<OneD, int> NumShape(2,0);
1139 
1140  for(int i = 0; i < GetExpSize(); ++i)
1141  {
1142  if((*m_exp)[i]->DetShapeType() == LibUtilities::eTriangle)
1143  {
1144  NumShape[0] += 1;
1145  }
1146  else // Quadrilateral element
1147  {
1148  NumShape[1] += 1;
1149  }
1150  }
1151 
1153  ::AllocateSharedPtr(m_session,2,NumShape);
1154  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:1052
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2046
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1036
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:969
void Nektar::MultiRegions::ExpList2D::v_SetUpPhysNormals ( )
privatevirtual

Set up the normals on each expansion.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1123 of file ExpList2D.cpp.

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

1124  {
1125  int i, j;
1126  for (i = 0; i < m_exp->size(); ++i)
1127  {
1128  for (j = 0; j < (*m_exp)[i]->GetNedges(); ++j)
1129  {
1130  (*m_exp)[i]->ComputeEdgeNormal(j);
1131  }
1132  }
1133  }
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1036
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 941 of file ExpList2D.cpp.

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

946  {
947  int i,j,f_npoints,offset;
948 
949  // Process each expansion.
950  for(i = 0; i < m_exp->size(); ++i)
951  {
952  // Get the number of points and the data offset.
953  f_npoints = (*m_exp)[i]->GetNumPoints(0)*
954  (*m_exp)[i]->GetNumPoints(1);
955  offset = m_phys_offset[i];
956 
957  // Process each point in the expansion.
958  for(j = 0; j < f_npoints; ++j)
959  {
960  // Upwind based on one-dimensional velocity.
961  if(Vn[offset + j] > 0.0)
962  {
963  Upwind[offset + j] = Fwd[offset + j];
964  }
965  else
966  {
967  Upwind[offset + j] = Bwd[offset + j];
968  }
969  }
970  }
971  }
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1036
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:1050
void Nektar::MultiRegions::ExpList2D::v_WriteVtkPieceHeader ( std::ostream &  outfile,
int  expansion,
int  istrip 
)
privatevirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1156 of file ExpList2D.cpp.

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