Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes | List of all members
Nektar::MultiRegions::ExpList1D Class Reference

This class is the abstraction of a one-dimensional multi-elemental expansions which is merely a collection of local expansions. More...

#include <ExpList1D.h>

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

Public Member Functions

 ExpList1D ()
 The default constructor. More...
 
 ExpList1D (const ExpList1D &In, const bool DeclareCoeffPhysArrays=true)
 The copy constructor. More...
 
 ExpList1D (const ExpList1D &In, const std::vector< unsigned int > &eIDs, const bool DeclareCoeffPhysArrays=true)
 Constructor copying only elements defined in eIds. More...
 
 ExpList1D (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::BasisKey &Ba, const SpatialDomains::MeshGraphSharedPtr &graph1D)
 Construct an ExpList1D from a given graph. More...
 
 ExpList1D (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &graph1D, const bool DeclareCoeffPhysArrays=true)
 This constructor sets up a list of local expansions based on an input graph1D. More...
 
 ExpList1D (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &graph1D, const SpatialDomains::CompositeMap &domain, const bool DeclareCoeffPhysArrays=true, const std::string var="DefaultVar", bool SetToOneSpaceDimension=false)
 This constructor sets up a list of local expansions based on an input compositeMap. More...
 
 ExpList1D (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::CompositeMap &domain, const SpatialDomains::MeshGraphSharedPtr &graph2D, const bool DeclareCoeffPhysArrays=true, const std::string variable="DefaultVar")
 Specialised constructor for Neumann boundary conditions in DisContField2D and ContField2D. More...
 
 ExpList1D (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::CompositeMap &domain, const SpatialDomains::MeshGraphSharedPtr &graph1D, int i, const bool DeclareCoeffPhysArrays=true)
 
 ExpList1D (const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, const ExpListSharedPtr > &bndConstraint, const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > &bndCond, const LocalRegions::ExpansionVector &locexp, const SpatialDomains::MeshGraphSharedPtr &graph2D, const PeriodicMap &periodicEdges, const bool DeclareCoeffPhysArrays=true, const std::string variable="DefaultVar")
 Specialised constructor for trace expansions. More...
 
virtual ~ExpList1D ()
 Destructor. More...
 
void PostProcess (LibUtilities::KernelSharedPtr kernel, Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray, NekDouble h, int elmId=0)
 Performs the post-processing on a specified element. More...
 
void PeriodicEval (Array< OneD, NekDouble > &inarray1, Array< OneD, NekDouble > &inarray2, NekDouble h, int nmodes, Array< OneD, NekDouble > &outarray)
 Evaluates the global spectral/hp expansion at some arbitray set of points. More...
 
void ParNormalSign (Array< OneD, NekDouble > &normsign)
 Set up the normals on each expansion. 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)
 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 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 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::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...
 
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)
 
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)
 
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 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 ()
 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 Array< OneD, NekDouble > > &Vec, 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 velocity field given by Vec. More...
 
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)
 Populate normals with the normals of all expansions. 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 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)
 
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_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)
 
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_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_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)
 
virtual void v_ClearGlobalLinSysManager (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. Sets up the storage for m_coeff and m_phys. More...
 
virtual void v_ReadGlobalOptimizationParameters ()
 
virtual void v_SetUpPhysNormals ()
 Set up the normals on each expansion. More...
 
virtual void v_WriteVtkPieceHeader (std::ostream &outfile, int expansion, int istrip)
 const StdRegions::StdExpansionVector &locexp); More...
 

Private Attributes

int m_firstIntEl
 
Array< OneD, NekDoublem_normSign
 

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...
 
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
 
boost::unordered_map< int, int > m_elmtToExpId
 Mapping from geometry ID of element to index inside m_exp. More...
 

Detailed Description

This class is the abstraction of a one-dimensional multi-elemental expansions 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}(x_i)=\sum_{e=1}^{{N_{\mathrm{el}}}} \sum_{n=0}^{N^{e}_m-1}\hat{u}_n^e\phi_n^e(x_i).\]

where ${N_{\mathrm{el}}}$ is the number of elements and $N^{e}_m$ is the local elemental number of expansion modes.

Instances of this class may be optionally constructed which use generalised segment expansions (LocalRegions::GenSegExp), rather than the standard segment expansions (LocalRegions::SegExp). LocalRegions::GenSegExp provides additional spatial data including segment normals and is enabled using the UseGenSegExp flag.

This class inherits all its variables and member functions from the base class MultiRegions::ExpList.

Definition at line 61 of file ExpList1D.h.

Constructor & Destructor Documentation

Nektar::MultiRegions::ExpList1D::ExpList1D ( )

The default constructor.

Assumes use of standard segment expansions only. All data storage areas are initialised to empty arrays by the default ExpList constructor.

Definition at line 77 of file ExpList1D.cpp.

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

77  :
78  ExpList()
79  {
80  SetExpType(e1D);
81  }
ExpList()
The default constructor.
Definition: ExpList.cpp:95
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:253
Nektar::MultiRegions::ExpList1D::ExpList1D ( const ExpList1D In,
const bool  DeclareCoeffPhysArrays = true 
)

The copy constructor.

Creates an identical copy of another ExpList1D object.

Definition at line 87 of file ExpList1D.cpp.

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

87  :
88  ExpList(In,DeclareCoeffPhysArrays)
89  {
90  SetExpType(e1D);
91  }
ExpList()
The default constructor.
Definition: ExpList.cpp:95
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:253
Nektar::MultiRegions::ExpList1D::ExpList1D ( const ExpList1D In,
const std::vector< unsigned int > &  eIDs,
const bool  DeclareCoeffPhysArrays = true 
)

Constructor copying only elements defined in eIds.

Definition at line 96 of file ExpList1D.cpp.

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

98  :
99  ExpList(In, eIDs, DeclareCoeffPhysArrays)
100  {
101  SetExpType(e1D);
102 
103  // Setup Default optimisation information.
104  int nel = GetExpSize();
107 
108  // Allocate storage for data and populate element offset lists.
110 
113  }
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:1001
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1899
void SetCoeffPhysOffsets(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
Definition: ExpList1D.cpp:709
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:2956
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:253
Nektar::MultiRegions::ExpList1D::ExpList1D ( const LibUtilities::SessionReaderSharedPtr pSession,
const LibUtilities::BasisKey Ba,
const SpatialDomains::MeshGraphSharedPtr graph1D 
)

Construct an ExpList1D from a given graph.

After initialising the data inherited through MultiRegions::ExpList, populate the expansion list from the segments defined in the supplied SpatialDomains::MeshGraph1D. All expansions in the graph are defined using the same LibUtilities::BasisKey which overrides that specified in graph1D.

See also
ExpList1D::ExpList1D(SpatialDomains::MeshGraph1D&, bool) for details.
Deprecated:
This constructor is no longer required since the basis key is now available from the graph.
Parameters
BaBasisKey describing quadrature points and number of modes.
graph1DDomain and expansion definitions.

Definition at line 131 of file ExpList1D.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::MultiRegions::ExpList::CreateCollections(), Nektar::MultiRegions::e1D, Nektar::MultiRegions::ExpList::GetExpSize(), 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().

133  :
134  ExpList(pSession,graph1D)
135  {
136  SetExpType(e1D);
137 
138  int id=0;
141 
142  const SpatialDomains::ExpansionMap &expansions
143  = graph1D->GetExpansions();
144 
145  // For each element in the mesh, create a segment expansion using
146  // the supplied BasisKey and segment geometry.
147  SpatialDomains::ExpansionMap::const_iterator expIt;
148  for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
149  {
150  if ((SegmentGeom = boost::dynamic_pointer_cast<
151  SpatialDomains::SegGeom>(
152  expIt->second->m_geomShPtr)))
153  {
155  ::AllocateSharedPtr(Ba,SegmentGeom);
156  seg->SetElmtId(id++);
157  (*m_exp).push_back(seg);
158  }
159  else
160  {
161  ASSERTL0(false,"dynamic cast to a SegGeom failed");
162  }
163  }
164 
165  // Setup Default optimisation information.
166  int nel = GetExpSize();
169 
170  // Allocate storage for data and populate element offset lists.
172 
173  m_coeffs = Array<OneD, NekDouble>(m_ncoeffs);
174  m_phys = Array<OneD, NekDouble>(m_npoints);
175 
178  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
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:1001
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:956
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:939
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1899
boost::shared_ptr< SegExp > SegExpSharedPtr
Definition: SegExp.h:266
boost::shared_ptr< SegGeom > SegGeomSharedPtr
Definition: Geometry2D.h:60
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:917
void SetCoeffPhysOffsets(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
Definition: ExpList1D.cpp:709
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:2956
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:253
std::map< int, ExpansionShPtr > ExpansionMap
Definition: MeshGraph.h:174
Nektar::MultiRegions::ExpList1D::ExpList1D ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr graph1D,
const bool  DeclareCoeffPhysArrays = true 
)

This constructor sets up a list of local expansions based on an input graph1D.

Given a mesh graph1D, 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 $x_i$ and local expansion coefficients $\hat{u}^e_n$ and allocates memory for the arrays m_coeffs and m_phys.

For each element its corresponding LibUtilities::BasisKey is retrieved and this is used to construct either a standard segment (LocalRegions::SegExp) or a generalised segment (LocalRegions::GenSegExp) which is stored in the list m_exp. Finally, ExpList::SetCoeffPhys is called to initialise the data storage areas and set up the offset arrays.

Parameters
graph1DA mesh, containing information about the domain and the spectral/hp element expansion.
UseGenSegExpIf true, create general segment expansions instead of just normal segment expansions.

Definition at line 201 of file ExpList1D.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::MultiRegions::ExpList::CreateCollections(), Nektar::MultiRegions::e1D, Nektar::MultiRegions::ExpList::GetExpSize(), 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().

203  :
204  ExpList(pSession,graph1D)
205  {
206  SetExpType(e1D);
207 
208  int id=0;
211 
212  // Retrieve the list of expansions
213  const SpatialDomains::ExpansionMap &expansions
214  = graph1D->GetExpansions();
215 
216  // Process each expansion in the graph
217  SpatialDomains::ExpansionMap::const_iterator expIt;
218  for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
219  {
220  // Retrieve basis key from expansion
221  LibUtilities::BasisKey bkey = expIt->second->m_basisKeyVector[0];
222 
223  if ((SegmentGeom = boost::dynamic_pointer_cast<
224  SpatialDomains::SegGeom>(
225  expIt->second->m_geomShPtr)))
226  {
228  ::AllocateSharedPtr(bkey, SegmentGeom);
229 
230  // Assign next ID
231  seg->SetElmtId(id++);
232 
233  // Add the expansion
234  (*m_exp).push_back(seg);
235  }
236  else
237  {
238  ASSERTL0(false,"dynamic cast to a SegGeom failed");
239  }
240  }
241 
242  // Setup Default optimisation information.
243  int nel = GetExpSize();
246 
247  // set up offset arrays.
249 
250  if(DeclareCoeffPhysArrays)
251  {
252  // Set up m_coeffs, m_phys.
253  m_coeffs = Array<OneD, NekDouble>(m_ncoeffs);
254  m_phys = Array<OneD, NekDouble>(m_npoints);
255  }
256 
259  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
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:1001
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:956
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:939
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1899
boost::shared_ptr< SegExp > SegExpSharedPtr
Definition: SegExp.h:266
boost::shared_ptr< SegGeom > SegGeomSharedPtr
Definition: Geometry2D.h:60
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:917
void SetCoeffPhysOffsets(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
Definition: ExpList1D.cpp:709
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:2956
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:253
std::map< int, ExpansionShPtr > ExpansionMap
Definition: MeshGraph.h:174
Nektar::MultiRegions::ExpList1D::ExpList1D ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr graph1D,
const SpatialDomains::CompositeMap domain,
const bool  DeclareCoeffPhysArrays = true,
const std::string  var = "DefaultVar",
bool  SetToOneSpaceDimension = false 
)

This constructor sets up a list of local expansions based on an input compositeMap.

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

For each element its corresponding LibUtilities::BasisKey is retrieved and this is used to construct either a standard segment (LocalRegions::SegExp) or a generalised segment (LocalRegions::GenSegExp) which is stored in the list m_exp. Finally, ExpList::SetCoeffPhys is called to initialise the data storage areas and set up the offset arrays.

Parameters
graph1DA mesh, containing information about the domain and the spectral/hp element expansion.
UseGenSegExpIf true, create general segment expansions instead of just normal segment expansions.

Definition at line 284 of file ExpList1D.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::MultiRegions::ExpList::CreateCollections(), Nektar::MultiRegions::ExpList::GetExpSize(), 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(), and SetCoeffPhysOffsets().

289  :
290  ExpList(pSession,graph1D)
291  {
292  int j,id=0;
296  SpatialDomains::CompositeMap::const_iterator compIt;
297 
298  // Retrieve the list of expansions
299  const SpatialDomains::ExpansionMap &expansions
300  = graph1D->GetExpansions(var);
301 
302  // Process each composite region in domain
303  for(compIt = domain.begin(); compIt != domain.end(); ++compIt)
304  {
305  comp = compIt->second;
306 
307  // Process each expansion in the graph
308  for(j = 0; j < compIt->second->size(); ++j)
309  {
310  SpatialDomains::ExpansionMap::const_iterator expIt;
311 
312  if((SegmentGeom = boost::dynamic_pointer_cast<
313  SpatialDomains::SegGeom>(
314  (*compIt->second)[j])))
315  {
316  // Retrieve basis key from expansion and define expansion
317  if((expIt = expansions.find(SegmentGeom->GetGlobalID())) != expansions.end())
318  {
319  LibUtilities::BasisKey bkey = expIt->second->m_basisKeyVector[0];
320 
321  if(SetToOneSpaceDimension)
322  {
323  SpatialDomains::SegGeomSharedPtr OneDSegmentGeom =
324  SegmentGeom->GenerateOneSpaceDimGeom();
325 
327  ::AllocateSharedPtr(bkey, OneDSegmentGeom);
328  }
329  else
330  {
332  ::AllocateSharedPtr(bkey, SegmentGeom);
333  }
334  }
335  else
336  {
337  ASSERTL0(false,"Failed to find basis key");
338  }
339 
340  }
341  else
342  {
343  ASSERTL0(false,"Failed to dynamic cast geometry to SegGeom");
344  }
345 
346 
347  // Assign next ID
348  seg->SetElmtId(id++);
349 
350  // Add the expansion
351  (*m_exp).push_back(seg);
352  }
353  }
354 
355  // Setup Default optimisation information.
356  int nel = GetExpSize();
359 
360  // set up offset arrays.
362 
363  if(DeclareCoeffPhysArrays)
364  {
365  // Set up m_coeffs, m_phys.
366  m_coeffs = Array<OneD, NekDouble>(m_ncoeffs);
367  m_phys = Array<OneD, NekDouble>(m_npoints);
368  }
369 
372  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
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:1001
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:956
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:939
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1899
boost::shared_ptr< SegExp > SegExpSharedPtr
Definition: SegExp.h:266
boost::shared_ptr< SegGeom > SegGeomSharedPtr
Definition: Geometry2D.h:60
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:917
boost::shared_ptr< GeometryVector > Composite
Definition: MeshGraph.h:114
void SetCoeffPhysOffsets(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
Definition: ExpList1D.cpp:709
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:2956
std::map< int, ExpansionShPtr > ExpansionMap
Definition: MeshGraph.h:174
Nektar::MultiRegions::ExpList1D::ExpList1D ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::CompositeMap domain,
const SpatialDomains::MeshGraphSharedPtr graph2D,
const bool  DeclareCoeffPhysArrays = true,
const std::string  variable = "DefaultVar" 
)

Specialised constructor for Neumann boundary conditions in DisContField2D and ContField2D.

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

See also
ExpList1D::ExpList1D(SpatialDomains::MeshGraph1D&, bool) for details.
Parameters
domainA domain, comprising of one or more composite regions,
graph2DA mesh, containing information about the domain and the spectral/hp element expansion.
UseGenSegExpIf true, create general segment expansions instead of just normal segment expansions.

Definition at line 388 of file ExpList1D.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::MultiRegions::ExpList::CreateCollections(), Nektar::MultiRegions::e1D, Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::MultiRegions::ExpList::m_coeffs, Nektar::MultiRegions::ExpList::m_globalOptParam, Nektar::MultiRegions::ExpList::m_graph, Nektar::MultiRegions::ExpList::m_ncoeffs, Nektar::MultiRegions::ExpList::m_npoints, Nektar::MultiRegions::ExpList::m_phys, SetCoeffPhysOffsets(), and Nektar::MultiRegions::ExpList::SetExpType().

392  :
393  ExpList(pSession,graph2D)
394  {
395  SetExpType(e1D);
396 
397  m_graph = graph2D;
398 
399  int j, id=0;
401  SpatialDomains::CompositeMap::const_iterator compIt;
404 
405  // Process each composite region.
406  for(compIt = domain.begin(); compIt != domain.end(); ++compIt)
407  {
408  comp = compIt->second;
409  // Process each expansion in the region.
410  for(j = 0; j < compIt->second->size(); ++j)
411  {
412  if((SegmentGeom = boost::dynamic_pointer_cast<
413  SpatialDomains::SegGeom>(
414  (*compIt->second)[j])))
415  {
416  // Retrieve the basis key from the expansion.
417  LibUtilities::BasisKey bkey
418  = boost::dynamic_pointer_cast<SpatialDomains::MeshGraph2D>(graph2D)->GetEdgeBasisKey(SegmentGeom, variable);
419 
421  ::AllocateSharedPtr(bkey, SegmentGeom);
422 
423  // Add the segment to the expansion list.
424  seg->SetElmtId(id++);
425  (*m_exp).push_back(seg);
426  }
427  else
428  {
429  ASSERTL0(false,"dynamic cast to a SegGeom failed");
430  }
431  }
432 
433  }
434 
435  // Setup Default optimisation information.
436  int nel = GetExpSize();
439 
440  // Allocate storage for data and populate element offset lists.
442 
443  // Set up m_coeffs, m_phys.
444  if(DeclareCoeffPhysArrays)
445  {
446  m_coeffs = Array<OneD, NekDouble>(m_ncoeffs);
447  m_phys = Array<OneD, NekDouble>(m_npoints);
448  }
449 
451  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
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:1001
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:956
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:939
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1899
boost::shared_ptr< SegExp > SegExpSharedPtr
Definition: SegExp.h:266
boost::shared_ptr< SegGeom > SegGeomSharedPtr
Definition: Geometry2D.h:60
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:917
SpatialDomains::MeshGraphSharedPtr m_graph
Mesh associated with this expansion list.
Definition: ExpList.h:913
boost::shared_ptr< GeometryVector > Composite
Definition: MeshGraph.h:114
void SetCoeffPhysOffsets(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
Definition: ExpList1D.cpp:709
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:2956
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:253
Nektar::MultiRegions::ExpList1D::ExpList1D ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::CompositeMap domain,
const SpatialDomains::MeshGraphSharedPtr graph1D,
int  i,
const bool  DeclareCoeffPhysArrays = true 
)
Nektar::MultiRegions::ExpList1D::ExpList1D ( const LibUtilities::SessionReaderSharedPtr pSession,
const Array< OneD, const ExpListSharedPtr > &  bndConstraint,
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > &  bndCond,
const LocalRegions::ExpansionVector locexp,
const SpatialDomains::MeshGraphSharedPtr graph2D,
const PeriodicMap periodicEdges,
const bool  DeclareCoeffPhysArrays = true,
const std::string  variable = "DefaultVar" 
)

Specialised constructor for trace expansions.

Store expansions for the trace space expansions used in DisContField2D.

Parameters
bndConstraintArray of ExpList1D objects each containing a 1D 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.
graph2D2D mesh corresponding to the expansion list.
periodicEdgesList of periodic edges.
UseGenSegExpIf true, create general segment expansions instead of just normal segment expansions.

Definition at line 469 of file ExpList1D.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::MultiRegions::ExpList::CreateCollections(), Nektar::MultiRegions::e1D, Nektar::SpatialDomains::eDirichlet, Nektar::LibUtilities::BasisKey::GetBasisType(), Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::LocalRegions::Expansion2D::GetGeom2D(), Nektar::StdRegions::StdExpansion::GetNedges(), 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().

477  :
478  ExpList(pSession,graph2D)
479  {
480  int i, j, id, elmtid = 0;
481  set<int> edgesDone;
482 
489 
490  SetExpType(e1D);
491 
492  map<int,int> EdgeDone;
493  map<int,int> NormalSet;
494 
496 
497  // First loop over boundary conditions to renumber
498  // Dirichlet boundaries
499  for(i = 0; i < bndCond.num_elements(); ++i)
500  {
501  if(bndCond[i]->GetBoundaryConditionType()
503  {
504  for(j = 0; j < bndConstraint[i]->GetExpSize(); ++j)
505  {
506  LibUtilities::BasisKey bkey = bndConstraint[i]
507  ->GetExp(j)->GetBasis(0)->GetBasisKey();
508  exp1D = bndConstraint[i]->GetExp(j)->
509  as<LocalRegions::Expansion1D>();
510  segGeom = exp1D->GetGeom1D();
511 
513  ::AllocateSharedPtr(bkey, segGeom);
514  edgesDone.insert(segGeom->GetEid());
515 
516  seg->SetElmtId(elmtid++);
517  (*m_exp).push_back(seg);
518  }
519  }
520  }
521 
523  LibUtilities::BasisKey> > edgeOrders;
525  LibUtilities::BasisKey> >::iterator it;
526 
527  for(i = 0; i < locexp.size(); ++i)
528  {
529  exp2D = locexp[i]->as<LocalRegions::Expansion2D>();
530 
531  for(j = 0; j < locexp[i]->GetNedges(); ++j)
532  {
533  segGeom = exp2D->GetGeom2D()->GetEdge(j);
534  id = segGeom->GetEid();
535  // Ignore Dirichlet edges
536  if (edgesDone.count(id) != 0)
537  {
538  continue;
539  }
540 
541  it = edgeOrders.find(id);
542 
543  if (it == edgeOrders.end())
544  {
545  edgeOrders.insert(std::make_pair(id, std::make_pair(
546  segGeom, locexp[i]->DetEdgeBasisKey(j))));
547  }
548  else // variable modes/points
549  {
550  LibUtilities::BasisKey edge
551  = locexp[i]->DetEdgeBasisKey(j);
552  LibUtilities::BasisKey existing
553  = it->second.second;
554 
555  int np1 = edge .GetNumPoints();
556  int np2 = existing.GetNumPoints();
557  int nm1 = edge .GetNumModes ();
558  int nm2 = existing.GetNumModes ();
559 
560  if (np2 >= np1 && nm2 >= nm1)
561  {
562  continue;
563  }
564  else if (np2 < np1 && nm2 < nm1)
565  {
566  it->second.second = edge;
567  }
568  else
569  {
570  ASSERTL0(false,
571  "inappropriate number of points/modes (max "
572  "num of points is not set with max order)");
573  }
574  }
575  }
576  }
577 
579  pSession->GetComm()->GetRowComm();
580  int nproc = vComm->GetSize(); // number of processors
581  int edgepr = vComm->GetRank(); // ID processor
582 
583  if (nproc > 1)
584  {
585  int eCnt = 0;
586 
587  // Count the number of edges on each partition
588  for(i = 0; i < locexp.size(); ++i)
589  {
590  eCnt += locexp[i]->GetNedges();
591  }
592 
593  // Set up the offset and the array that will contain the list of
594  // edge IDs, then reduce this across processors.
595  Array<OneD, int> edgesCnt(nproc, 0);
596  edgesCnt[edgepr] = eCnt;
597  vComm->AllReduce(edgesCnt, LibUtilities::ReduceSum);
598 
599  // Set up offset array.
600  int totEdgeCnt = Vmath::Vsum(nproc, edgesCnt, 1);
601  Array<OneD, int> eTotOffsets(nproc,0);
602  for (i = 1; i < nproc; ++i)
603  {
604  eTotOffsets[i] = eTotOffsets[i-1] + edgesCnt[i-1];
605  }
606 
607  // Local list of the edges per element
608  Array<OneD, int> EdgesTotID(totEdgeCnt, 0);
609  Array<OneD, int> EdgesTotNm(totEdgeCnt, 0);
610  Array<OneD, int> EdgesTotPnts(totEdgeCnt, 0);
611 
612  int cntr = eTotOffsets[edgepr];
613 
614  for(i = 0; i < locexp.size(); ++i)
615  {
616  exp2D = locexp[i]->as<LocalRegions::Expansion2D>();
617 
618  int nedges = locexp[i]->GetNedges();
619 
620  for(j = 0; j < nedges; ++j, ++cntr)
621  {
622  LibUtilities::BasisKey bkeyEdge =
623  locexp[i]->DetEdgeBasisKey(j);
624  EdgesTotID [cntr] = exp2D->GetGeom2D()->GetEid(j);
625  EdgesTotNm [cntr] = bkeyEdge.GetNumModes();
626  EdgesTotPnts[cntr] = bkeyEdge.GetNumPoints();
627  }
628  }
629 
630  vComm->AllReduce(EdgesTotID, LibUtilities::ReduceSum);
631  vComm->AllReduce(EdgesTotNm, LibUtilities::ReduceSum);
632  vComm->AllReduce(EdgesTotPnts, LibUtilities::ReduceSum);
633 
634  for (i = 0; i < totEdgeCnt; ++i)
635  {
636  it = edgeOrders.find(EdgesTotID[i]);
637 
638  if (it == edgeOrders.end())
639  {
640  continue;
641  }
642 
643  LibUtilities::BasisKey existing
644  = it->second.second;
645  LibUtilities::BasisKey edge(
646  existing.GetBasisType(), EdgesTotNm[i],
647  LibUtilities::PointsKey(EdgesTotPnts[i],
648  existing.GetPointsType()));
649 
650 
651  int np1 = edge .GetNumPoints();
652  int np2 = existing.GetNumPoints();
653  int nm1 = edge .GetNumModes ();
654  int nm2 = existing.GetNumModes ();
655 
656  if (np2 >= np1 && nm2 >= nm1)
657  {
658  continue;
659  }
660  else if (np2 < np1 && nm2 < nm1)
661  {
662  it->second.second = edge;
663  }
664  else
665  {
666  ASSERTL0(false,
667  "inappropriate number of points/modes (max "
668  "num of points is not set with max order)");
669  }
670  }
671  }
672 
673  for (it = edgeOrders.begin(); it != edgeOrders.end(); ++it)
674  {
676  ::AllocateSharedPtr(it->second.second, it->second.first);
677  seg->SetElmtId(elmtid++);
678  (*m_exp).push_back(seg);
679  }
680 
681  // Setup Default optimisation information.
682  int nel = GetExpSize();
685 
686  // Set up offset information and array sizes
688 
689  // Set up m_coeffs, m_phys.
690  if(DeclareCoeffPhysArrays)
691  {
692  m_coeffs = Array<OneD, NekDouble>(m_ncoeffs);
693  m_phys = Array<OneD, NekDouble>(m_npoints);
694  }
695 
697  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
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:1001
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:956
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:939
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1899
boost::shared_ptr< SegExp > SegExpSharedPtr
Definition: SegExp.h:266
boost::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:53
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:917
void SetCoeffPhysOffsets(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
Definition: ExpList1D.cpp:709
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< Geometry1D > Geometry1DSharedPtr
Definition: Geometry1D.h:48
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:2956
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:253
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:723
boost::shared_ptr< Expansion1D > Expansion1DSharedPtr
Definition: Expansion1D.h:53
boost::shared_ptr< Expansion2D > Expansion2DSharedPtr
Definition: Expansion1D.h:49
Nektar::MultiRegions::ExpList1D::~ExpList1D ( )
virtual

Destructor.

Definition at line 733 of file ExpList1D.cpp.

734  {
735  }

Member Function Documentation

void Nektar::MultiRegions::ExpList1D::ParNormalSign ( Array< OneD, NekDouble > &  normsign)

Set up the normals on each expansion.

void Nektar::MultiRegions::ExpList1D::PeriodicEval ( Array< OneD, NekDouble > &  inarray1,
Array< OneD, NekDouble > &  inarray2,
NekDouble  h,
int  nmodes,
Array< OneD, NekDouble > &  outarray 
)

Evaluates the global spectral/hp expansion at some arbitray set of points.

Given the elemental coefficients $\hat{u}_n^e$ of an expansion, periodically evaluate the spectral/hp expansion $u^{\delta}(\boldsymbol{x})$ at arbitrary points.

Parameters
inarray1An array of size $N_{\mathrm{eof}}$ containing the local coefficients $\hat{u}_n^e$.
inarray2Contains the set of evaluation points.
hThe mesh spacing.
nmodesThe number of polynomial modes for each element (we consider that each element has the same number of polynomial modes).
outarrayContains the resulting values at the evaluation points

Definition at line 846 of file ExpList1D.cpp.

References Nektar::MultiRegions::ExpList::GetExpSize(), and Polylib::jacobfd().

Referenced by PostProcess().

850  {
851  int i,j,r;
852 
853  // Get the number of elements in the domain
854  int num_elm = GetExpSize();
855 
856  // initializing the outarray
857  for(i = 0; i < outarray.num_elements(); i++)
858  {
859  outarray[i] = 0.0;
860  }
861 
862  // Make a copy for further modification
863  int x_size = inarray2.num_elements();
864  Array<OneD,NekDouble> x_values_cp(x_size);
865 
866  // Determining the element to which the x belongs
867  Array<OneD,int> x_elm(x_size);
868  for(i = 0; i < x_size; i++ )
869  {
870  x_elm[i] = (int)floor(inarray2[i]/h);
871  }
872 
873  // Clamp indices periodically
874  for(i = 0; i < x_size; i++)
875  {
876  while(x_elm[i] < 0)
877  {
878  x_elm[i] += num_elm;
879  }
880  while(x_elm[i] >= num_elm)
881  {
882  x_elm[i] -= num_elm ;
883  }
884  }
885 
886  // Map the values of x to [-1 1] on its interval
887  for(i = 0; i < x_size; i++)
888  {
889  x_values_cp[i] = (inarray2[i]/h - floor(inarray2[i]/h))*2 - 1.0;
890  }
891 
892  // Evaluate the jocobi polynomials
893  // (Evaluating the base at some points other than the quadrature
894  // points). Should it be added to the base class????
895  Array<TwoD,NekDouble> jacobi_poly(nmodes,x_size);
896  for(i = 0; i < nmodes; i++)
897  {
898  Polylib::jacobfd(x_size,x_values_cp.get(),
899  jacobi_poly.get()+i*x_size,NULL,i,0.0,0.0);
900  }
901 
902  // Evaluate the function values
903  for(r = 0; r < nmodes; r++)
904  {
905  for(j = 0; j < x_size; j++)
906  {
907  int index = ((x_elm[j])*nmodes)+r;
908  outarray[j] += inarray1[index]*jacobi_poly[r][j];
909  }
910  }
911 
912  }
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1899
void jacobfd(const int np, const double *z, double *poly_in, double *polyd, const int n, const double alpha, const double beta)
Routine to calculate Jacobi polynomials, , and their first derivative, .
Definition: Polylib.cpp:1946
void Nektar::MultiRegions::ExpList1D::PostProcess ( LibUtilities::KernelSharedPtr  kernel,
Array< OneD, NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
NekDouble  h,
int  elmId = 0 
)

Performs the post-processing on a specified element.

To perform post-processing on the entire domain use elmtId = 0.

Parameters
kernelThe post-processing kernel.
inarrayThe set of evaluation points.
outarrayContains the resulting post-processed solution for element elmId.
hThe mesh spacing.
elmIdOptionally specifies which element to perform the post-processing on (0=whole domain).

Definition at line 747 of file ExpList1D.cpp.

References Nektar::MultiRegions::ExpList::GetCoeffs(), Nektar::MultiRegions::ExpList::GetExp(), PeriodicEval(), and Nektar::LibUtilities::PointsManager().

753  {
754  int i,j,r;
755 
756  // get the local element expansion of the elmId element
758 
759  // Get the quadrature points and weights required for integration
760  int quad_npoints = elmExp->GetTotPoints();
761  LibUtilities::PointsKey quadPointsKey(quad_npoints,
762  elmExp->GetPointsType(0));
763  Array<OneD,NekDouble> quad_points
764  = LibUtilities::PointsManager()[quadPointsKey]->GetZ();
765  Array<OneD,NekDouble> quad_weights
766  = LibUtilities::PointsManager()[quadPointsKey]->GetW();
767 
768  // Declare variable for the local kernel breaks
769  int kernel_width = kernel->GetKernelWidth();
770  Array<OneD,NekDouble> local_kernel_breaks(kernel_width+1);
771 
772  // Declare variable for the transformed quadrature points
773  Array<OneD,NekDouble> mapped_quad_points(quad_npoints);
774 
775  // For each evaluation point
776  for(i = 0; i < inarray.num_elements(); i++)
777  {
778  // Move the center of the kernel to the current point
779  kernel->MoveKernelCenter(inarray[i],local_kernel_breaks);
780 
781  // Find the mesh breaks under the kernel support
782  Array<OneD,NekDouble> mesh_breaks;
783  kernel->FindMeshUnderKernel(local_kernel_breaks,h,mesh_breaks);
784 
785  // Sort the total breaks for integration purposes
786  int total_nbreaks = local_kernel_breaks.num_elements() +
787  mesh_breaks.num_elements();
788  // number of the total breaks
789  Array<OneD,NekDouble> total_breaks(total_nbreaks);
790  kernel->Sort(local_kernel_breaks,mesh_breaks,total_breaks);
791 
792  // Integrate the product of kernel and function over the total
793  // breaks
794  NekDouble integral_value = 0.0;
795  for(j = 0; j < total_breaks.num_elements()-1; j++)
796  {
797  NekDouble a = total_breaks[j];
798  NekDouble b = total_breaks[j+1];
799 
800  // Map the quadrature points to the appropriate interval
801  for(r = 0; r < quad_points.num_elements(); r++)
802  {
803  mapped_quad_points[r]
804  = (quad_points[r] + 1.0) * 0.5 * (b - a) + a;
805  }
806 
807  // Evaluate the function at the transformed quadrature
808  // points
809  Array<OneD,NekDouble> u_value(quad_npoints);
810  Array<OneD,NekDouble> coeffs = GetCoeffs();
811 
812  PeriodicEval(coeffs,mapped_quad_points,h,
813  elmExp->GetBasisNumModes(0),u_value);
814 
815  // Evaluate the kernel at the transformed quadrature points
816  Array<OneD,NekDouble> k_value(quad_npoints);
817  kernel->EvaluateKernel(mapped_quad_points,h,k_value);
818 
819  // Integrate
820  for(r = 0; r < quad_npoints; r++)
821  {
822  integral_value += (b - a) * 0.5 * k_value[r]
823  * u_value[r] * quad_weights[r];
824  }
825  }
826  outarray[i] = integral_value/h;
827  }
828  }
const Array< OneD, const NekDouble > & GetCoeffs() const
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expa...
Definition: ExpList.h:1837
const boost::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:1920
void PeriodicEval(Array< OneD, NekDouble > &inarray1, Array< OneD, NekDouble > &inarray2, NekDouble h, int nmodes, Array< OneD, NekDouble > &outarray)
Evaluates the global spectral/hp expansion at some arbitray set of points.
Definition: ExpList1D.cpp:846
PointsManagerT & PointsManager(void)
double NekDouble
boost::shared_ptr< StdExpansion > StdExpansionSharedPtr
void Nektar::MultiRegions::ExpList1D::SetCoeffPhysOffsets ( void  )
private

Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m_coeff and m_phys.

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 also 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 709 of file ExpList1D.cpp.

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

710  {
711  int i;
712 
713  // Set up offset information and array sizes
714  m_coeff_offset = Array<OneD,int>(m_exp->size());
715  m_phys_offset = Array<OneD,int>(m_exp->size());
716  m_offset_elmt_id = Array<OneD,int>(m_exp->size());
717 
718  m_ncoeffs = m_npoints = 0;
719 
720  for(i = 0; i < m_exp->size(); ++i)
721  {
723  m_phys_offset [i] = m_npoints;
724  m_offset_elmt_id[i] = i;
725  m_ncoeffs += (*m_exp)[i]->GetNcoeffs();
726  m_npoints += (*m_exp)[i]->GetTotPoints();
727  }
728  }
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:988
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:977
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:991
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:917
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:999
void Nektar::MultiRegions::ExpList1D::v_GetNormals ( Array< OneD, Array< OneD, NekDouble > > &  normals)
protectedvirtual

Populate normals with the normals of all expansions.

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 1074 of file ExpList1D.cpp.

References ASSERTL1, Nektar::MultiRegions::ExpList::GetCoordim(), Nektar::LocalRegions::Expansion1D::GetLeftAdjacentElementExp(), Nektar::LibUtilities::Interp1D(), Nektar::MultiRegions::ExpList::m_exp, and Nektar::MultiRegions::ExpList::m_phys_offset.

1076  {
1077  int i,j,k,e_npoints,offset;
1079  Array<OneD,Array<OneD,NekDouble> > locnormals;
1080  Array<OneD,Array<OneD,NekDouble> > locnormals2;
1081  Array<OneD,Array<OneD,NekDouble> > Norms;
1082  // Assume whole array is of same coordinate dimension
1083  int coordim = GetCoordim(0);
1084 
1085  ASSERTL1(normals.num_elements() >= coordim,
1086  "Output vector does not have sufficient dimensions to "
1087  "match coordim");
1088 
1089  for (i = 0; i < m_exp->size(); ++i)
1090  {
1091  LocalRegions::Expansion1DSharedPtr loc_exp = (*m_exp)[i]->as<LocalRegions::Expansion1D>();
1092 
1094  loc_exp->GetLeftAdjacentElementExp();
1095 
1096  int edgeNumber = loc_exp->GetLeftAdjacentElementEdge();
1097 
1098  // Get the number of points and normals for this expansion.
1099  e_npoints = (*m_exp)[i]->GetNumPoints(0);
1100 
1101  locnormals = loc_elmt->GetEdgeNormal(edgeNumber);
1102  int e_nmodes = loc_exp->GetBasis(0)->GetNumModes();
1103  int loc_nmodes = loc_elmt->GetBasis(0)->GetNumModes();
1104 
1105  if (e_nmodes != loc_nmodes)
1106  {
1107  if (loc_exp->GetRightAdjacentElementEdge() >= 0)
1108  {
1110  loc_exp->GetRightAdjacentElementExp();
1111 
1112  int EdgeNumber = loc_exp->GetRightAdjacentElementEdge();
1113  // Serial case: right element is connected so we can
1114  // just grab that normal.
1115  locnormals = loc_elmt->GetEdgeNormal(EdgeNumber);
1116 
1117  offset = m_phys_offset[i];
1118 
1119  // Process each point in the expansion.
1120  for (j = 0; j < e_npoints; ++j)
1121  {
1122  // Process each spatial dimension and copy the values
1123  // into the output array.
1124  for (k = 0; k < coordim; ++k)
1125  {
1126  normals[k][offset + j] = -locnormals[k][j];
1127  }
1128  }
1129  }
1130  else
1131  {
1132  // Parallel case: need to interpolate normal.
1133  Array<OneD, Array<OneD, NekDouble> > normal(coordim);
1134 
1135  for (int p = 0; p < coordim; ++p)
1136  {
1137  normal[p] = Array<OneD, NekDouble>(e_npoints,0.0);
1138  LibUtilities::PointsKey to_key =
1139  loc_exp->GetBasis(0)->GetPointsKey();
1140  LibUtilities::PointsKey from_key =
1141  loc_elmt->GetBasis(0)->GetPointsKey();
1142  LibUtilities::Interp1D(from_key,
1143  locnormals[p],
1144  to_key,
1145  normal[p]);
1146  }
1147 
1148  offset = m_phys_offset[i];
1149 
1150  // Process each point in the expansion.
1151  for (j = 0; j < e_npoints; ++j)
1152  {
1153  // Process each spatial dimension and copy the values
1154  // into the output array.
1155  for (k = 0; k < coordim; ++k)
1156  {
1157  normals[k][offset + j] = normal[k][j];
1158  }
1159  }
1160  }
1161  }
1162  else
1163  {
1164  // Get the physical data offset for this expansion.
1165  offset = m_phys_offset[i];
1166 
1167  // Process each point in the expansion.
1168  for (j = 0; j < e_npoints; ++j)
1169  {
1170  // Process each spatial dimension and copy the values
1171  // into the output array.
1172  for (k = 0; k < coordim; ++k)
1173  {
1174  normals[k][offset + j] = locnormals[k][j];
1175  }
1176  }
1177  }
1178  }
1179  }
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:977
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:991
void Interp1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
this function interpolates a 1D function evaluated at the quadrature points of the basis fbasis0 to ...
Definition: Interp.cpp:54
boost::shared_ptr< Geometry1D > Geometry1DSharedPtr
Definition: Geometry1D.h:48
int GetCoordim(int eid)
This function returns the dimension of the coordinates of the element eid.
Definition: ExpList.h:1797
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
boost::shared_ptr< Expansion1D > Expansion1DSharedPtr
Definition: Expansion1D.h:53
boost::shared_ptr< Expansion2D > Expansion2DSharedPtr
Definition: Expansion1D.h:49
void Nektar::MultiRegions::ExpList1D::v_ReadGlobalOptimizationParameters ( )
privatevirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1184 of file ExpList1D.cpp.

1185  {
1186 // Array<OneD, int> NumShape(1,0);
1187 // NumShape[0] = GetExpSize();
1188 //
1189 // int one = 1;
1190 // m_globalOptParam = MemoryManager<NekOptimize::GlobalOptParam>
1191 // ::AllocateSharedPtr(m_session,one,NumShape);
1192  }
void Nektar::MultiRegions::ExpList1D::v_SetUpPhysNormals ( )
privatevirtual

Set up the normals on each expansion.

Sets up the normals on all edges of expansions in the domain.

Parameters
locexpComplete list of domain expansions.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 949 of file ExpList1D.cpp.

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

950  {
951  int i, j;
952  for (i = 0; i < m_exp->size(); ++i)
953  {
954  for (j = 0; j < (*m_exp)[i]->GetNverts(); ++j)
955  {
956  (*m_exp)[i]->ComputeVertexNormal(j);
957  }
958  }
959  }
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:977
void Nektar::MultiRegions::ExpList1D::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 
)
protectedvirtual

Upwind the Fwd and Bwd states based on the velocity field given by Vec.

Upwind the left and right states given by the Arrays Fwd and Bwd using the vector quantity Vec and ouput the upwinded value in the array upwind.

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

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 971 of file ExpList1D.cpp.

References ASSERTL1, Nektar::MultiRegions::ExpList::GetCoordim(), Nektar::MultiRegions::ExpList::m_exp, and Nektar::MultiRegions::ExpList::m_phys_offset.

976  {
977  int i,j,k,e_npoints,offset;
978  Array<OneD,NekDouble> normals;
979  NekDouble Vn;
980 
981  // Assume whole array is of same coordimate dimension
982  int coordim = GetCoordim(0);
983 
984  ASSERTL1(Vec.num_elements() >= coordim,
985  "Input vector does not have sufficient dimensions to "
986  "match coordim");
987 
988  // Process each expansion
989  for(i = 0; i < m_exp->size(); ++i)
990  {
991  // Get the number of points in the expansion and the normals.
992  e_npoints = (*m_exp)[i]->GetNumPoints(0);
993  normals = (*m_exp)[i]->GetPhysNormals();
994 
995  // Get the physical data offset of the expansion in m_phys.
996  offset = m_phys_offset[i];
997 
998  // Compute each data point.
999  for(j = 0; j < e_npoints; ++j)
1000  {
1001  // Calculate normal velocity.
1002  Vn = 0.0;
1003  for(k = 0; k < coordim; ++k)
1004  {
1005  Vn += Vec[k][offset+j]*normals[k*e_npoints + j];
1006  }
1007 
1008  // Upwind based on direction of normal velocity.
1009  if(Vn > 0.0)
1010  {
1011  Upwind[offset + j] = Fwd[offset + j];
1012  }
1013  else
1014  {
1015  Upwind[offset + j] = Bwd[offset + j];
1016  }
1017  }
1018  }
1019  }
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:977
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:991
double NekDouble
int GetCoordim(int eid)
This function returns the dimension of the coordinates of the element eid.
Definition: ExpList.h:1797
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
void Nektar::MultiRegions::ExpList1D::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.

See also
ExpList1D::Upwind( const Array<OneD, const Array<OneD, NekDouble> >, const Array<OneD, const NekDouble>, const Array<OneD, const NekDouble>, Array<OneD, NekDouble>, int)
Parameters
VnVelocity field.
FwdLeft state.
BwdRight state.
UpwindOutput vector.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1034 of file ExpList1D.cpp.

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

1039  {
1040  int i,j,e_npoints,offset;
1041  Array<OneD,NekDouble> normals;
1042 
1043  // Process each expansion.
1044  for(i = 0; i < m_exp->size(); ++i)
1045  {
1046  // Get the number of points and the data offset.
1047  e_npoints = (*m_exp)[i]->GetNumPoints(0);
1048  offset = m_phys_offset[i];
1049 
1050  // Process each point in the expansion.
1051  for(j = 0; j < e_npoints; ++j)
1052  {
1053  // Upwind based on one-dimensional velocity.
1054  if(Vn[offset + j] > 0.0)
1055  {
1056  Upwind[offset + j] = Fwd[offset + j];
1057  }
1058  else
1059  {
1060  Upwind[offset + j] = Bwd[offset + j];
1061  }
1062  }
1063  }
1064  }
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:977
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:991
void Nektar::MultiRegions::ExpList1D::v_WriteVtkPieceHeader ( std::ostream &  outfile,
int  expansion,
int  istrip 
)
privatevirtual

const StdRegions::StdExpansionVector &locexp);

Writes out the header for a <PIECE> VTK XML segment describing the geometric information which comprises this element. This includes vertex coordinates for each quadrature point, vertex connectivity information, cell types and cell offset data.

Parameters
outfileOutput stream to write data to.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1203 of file ExpList1D.cpp.

1204  {
1205  int i,j;
1206  int nquad0 = (*m_exp)[expansion]->GetNumPoints(0);
1207  int ntot = nquad0;
1208  int ntotminus = (nquad0-1);
1209 
1210  Array<OneD,NekDouble> coords[3];
1211  coords[0] = Array<OneD,NekDouble>(ntot, 0.0);
1212  coords[1] = Array<OneD,NekDouble>(ntot, 0.0);
1213  coords[2] = Array<OneD,NekDouble>(ntot, 0.0);
1214  (*m_exp)[expansion]->GetCoords(coords[0],coords[1],coords[2]);
1215 
1216  outfile << " <Piece NumberOfPoints=\""
1217  << ntot << "\" NumberOfCells=\""
1218  << ntotminus << "\">" << endl;
1219  outfile << " <Points>" << endl;
1220  outfile << " <DataArray type=\"Float64\" "
1221  << "NumberOfComponents=\"3\" format=\"ascii\">" << endl;
1222  outfile << " ";
1223  for (i = 0; i < ntot; ++i)
1224  {
1225  for (j = 0; j < 3; ++j)
1226  {
1227  outfile << setprecision(8) << scientific
1228  << (float)coords[j][i] << " ";
1229  }
1230  outfile << endl;
1231  }
1232  outfile << endl;
1233  outfile << " </DataArray>" << endl;
1234  outfile << " </Points>" << endl;
1235  outfile << " <Cells>" << endl;
1236  outfile << " <DataArray type=\"Int32\" "
1237  << "Name=\"connectivity\" format=\"ascii\">" << endl;
1238  for (i = 0; i < nquad0-1; ++i)
1239  {
1240  outfile << i << " " << i+1 << endl;
1241  }
1242  outfile << endl;
1243  outfile << " </DataArray>" << endl;
1244  outfile << " <DataArray type=\"Int32\" "
1245  << "Name=\"offsets\" format=\"ascii\">" << endl;
1246  for (i = 0; i < ntotminus; ++i)
1247  {
1248  outfile << i*2+2 << " ";
1249  }
1250  outfile << endl;
1251  outfile << " </DataArray>" << endl;
1252  outfile << " <DataArray type=\"UInt8\" "
1253  << "Name=\"types\" format=\"ascii\">" << endl;
1254  for (i = 0; i < ntotminus; ++i)
1255  {
1256  outfile << "3 ";
1257  }
1258  outfile << endl;
1259  outfile << " </DataArray>" << endl;
1260  outfile << " </Cells>" << endl;
1261  outfile << " <PointData>" << endl;
1262  }

Member Data Documentation

int Nektar::MultiRegions::ExpList1D::m_firstIntEl
private

Definition at line 190 of file ExpList1D.h.

Array<OneD, NekDouble> Nektar::MultiRegions::ExpList1D::m_normSign
private

Definition at line 192 of file ExpList1D.h.