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.
 ExpList1D (const ExpList1D &In, const bool DeclareCoeffPhysArrays=true)
 The copy constructor.
 ExpList1D (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::BasisKey &Ba, const SpatialDomains::MeshGraphSharedPtr &graph1D)
 Construct an ExpList1D from a given graph.
 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.
 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.
 ExpList1D (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.
 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.
virtual ~ExpList1D ()
 Destructor.
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.
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.
void ParNormalSign (Array< OneD, NekDouble > &normsign)
 Set up the normals on each expansion.
- Public Member Functions inherited from Nektar::MultiRegions::ExpList
 ExpList ()
 The default constructor.
 ExpList (const LibUtilities::SessionReaderSharedPtr &pSession)
 The default constructor.
 ExpList (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 The default constructor.
 ExpList (const ExpList &in, const bool DeclareCoeffPhysArrays=true)
 The copy constructor.
virtual ~ExpList ()
 The default destructor.
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$.
int GetNcoeffs (const int eid) const
 Returns the total number of local degrees of freedom for element eid.
ExpansionType GetExpType (void)
 Returns the type of the expansion.
void SetExpType (ExpansionType Type)
 Returns the type of the expansion.
int EvalBasisNumModesMax (void) const
 Evaulates the maximum number of modes in the elemental basis order over all elements.
const Array< OneD, int > EvalBasisNumModesMaxPerExp (void) const
 Returns the vector of the number of modes in the elemental basis order over all elements.
int GetTotPoints (void) const
 Returns the total number of quadrature points m_npoints $=Q_{\mathrm{tot}}$.
int GetTotPoints (const int eid) const
 Returns the total number of quadrature points for eid's element $=Q_{\mathrm{tot}}$.
int GetNpoints (void) const
 Returns the total number of quadrature points m_npoints $=Q_{\mathrm{tot}}$.
int Get1DScaledTotPoints (const NekDouble scale) const
 Returns the total number of qudature points scaled by the factor scale on each 1D direction.
void SetWaveSpace (const bool wavespace)
 Sets the wave space to the one of the possible configuration true or false.
void SetModifiedBasis (const bool modbasis)
 Set Modified Basis for the stability analysis.
void SetPhys (int i, NekDouble val)
 Set the i th value of m_phys to value val.
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.
void SetPhys (const Array< OneD, const NekDouble > &inarray)
 Fills the array m_phys.
void SetPhysArray (Array< OneD, NekDouble > &inarray)
 Sets the array m_phys.
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.
bool GetPhysState (void) const
 This function indicates whether the array of physical values $\boldsymbol{u}_l$ (implemented as m_phys) is filled or not.
NekDouble PhysIntegral (void)
 This function integrates a function $f(\boldsymbol{x})$ over the domain consisting of all the elements of the expansion.
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.
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})$.
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.
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.
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.
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.
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.
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.
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.
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.
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$.
void HomogeneousFwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
void HomogeneousBwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
void DealiasedProd (const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
void GetBCValues (Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
void NormVectorIProductWRTBase (Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
void ApplyGeomInfo ()
 Apply geometry information to each expansion.
void WriteTecplotHeader (std::ofstream &outfile, std::string var="")
void WriteTecplotZone (std::ofstream &outfile, int expansion=-1)
void WriteTecplotField (std::ofstream &outfile, int expansion=-1)
void WriteTecplotConnectivity (std::ofstream &outfile, int expansion=-1)
void WriteVtkHeader (std::ofstream &outfile)
void WriteVtkFooter (std::ofstream &outfile)
void WriteVtkPieceHeader (std::ofstream &outfile, int expansion)
void WriteVtkPieceFooter (std::ofstream &outfile, int expansion)
void WriteVtkPieceData (std::ofstream &outfile, int expansion, std::string var="v")
int GetCoordim (int eid)
 This function returns the dimension of the coordinates of the element eid.
void SetCoeff (int i, NekDouble val)
 Set the i th coefficiient in m_coeffs to value val.
void SetCoeffs (int i, NekDouble val)
 Set the i th coefficiient in m_coeffs to value val.
void SetCoeffsArray (Array< OneD, NekDouble > &inarray)
 Set the m_coeffs array to inarray.
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.
void ImposeDirichletConditions (Array< OneD, NekDouble > &outarray)
 Impose Dirichlet Boundary Conditions onto Array.
void FillBndCondFromField (void)
 Fill Bnd Condition expansion from the values stored in expansion.
void LocalToGlobal (void)
 Put the coefficients into global ordering using m_coeffs.
void GlobalToLocal (void)
 Put the coefficients into local ordering and place in m_coeffs.
NekDouble GetCoeff (int i)
 Get the i th value (coefficient) of m_coeffs.
NekDouble GetCoeffs (int i)
 Get the i th value (coefficient) of m_coeffs.
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.
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.
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.
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.
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.
void SetHomo1DSpecVanVisc (Array< OneD, NekDouble > visc)
 This function sets the Spectral Vanishing Viscosity in homogeneous1D expansion.
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.
LibUtilities::TranspositionSharedPtr GetTransposition (void)
 This function returns the transposition class associaed with the homogeneous expansion.
NekDouble GetHomoLen (void)
 This function returns the Width of homogeneous direction associaed with the homogeneous expansion.
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.
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.
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.
int GetExpSize (void)
 This function returns the number of elements in the expansion.
int GetNumElmts (void)
 This function returns the number of elements in the expansion which may be different for a homogeoenous extended expansionp.
const boost::shared_ptr
< LocalRegions::ExpansionVector
GetExp () const
 This function returns the vector of elements in the expansion.
LocalRegions::ExpansionSharedPtrGetExp (int n) const
 This function returns (a shared pointer to) the local elemental expansion of the $n^{\mathrm{th}}$ element.
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.
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.
int GetPhys_Offset (int n) const
 Get the start offset position for a global list of m_phys correspoinding to element n.
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.
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.
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.
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.
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
const Array< OneD, const
boost::shared_ptr< ExpList > > & 
GetBndCondExpansions ()
boost::shared_ptr< ExpList > & UpdateBndCondExpansion (int i)
void Upwind (const Array< OneD, const Array< OneD, NekDouble > > &Vec, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
void Upwind (const Array< OneD, const NekDouble > &Vn, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
boost::shared_ptr< ExpList > & GetTrace ()
boost::shared_ptr
< AssemblyMapDG > & 
GetTraceMap (void)
const Array< OneD, const int > & GetTraceBndMap (void)
void GetNormals (Array< OneD, Array< OneD, NekDouble > > &normals)
void AddTraceIntegral (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
void AddTraceIntegral (const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
void AddFwdBwdTraceIntegral (const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
void GetFwdBwdTracePhys (Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
void GetFwdBwdTracePhys (const Array< OneD, const NekDouble > &field, Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
void ExtractTracePhys (Array< OneD, NekDouble > &outarray)
void ExtractTracePhys (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
const Array< OneD, const
SpatialDomains::BoundaryConditionShPtr > & 
GetBndConditions ()
Array< OneD,
SpatialDomains::BoundaryConditionShPtr > & 
UpdateBndConditions ()
void EvaluateBoundaryConditions (const NekDouble time=0.0, const std::string varName="", const NekDouble=NekConstants::kNekUnsetDouble, const NekDouble=NekConstants::kNekUnsetDouble)
void GeneralMatrixOp (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
 This function calculates the result of the multiplication of a matrix of type specified by mkey with a vector given by inarray.
void GeneralMatrixOp_IterPerExp (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void SetUpPhysNormals ()
void GetBoundaryToElmtMap (Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
void GeneralGetFieldDefinitions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef, int NumHomoDir=0, Array< OneD, LibUtilities::BasisSharedPtr > &HomoBasis=LibUtilities::NullBasisSharedPtr1DArray, std::vector< NekDouble > &HomoLen=LibUtilities::NullNekDoubleVector, std::vector< unsigned int > &HomoZIDs=LibUtilities::NullUnsignedIntVector, std::vector< unsigned int > &HomoYIDs=LibUtilities::NullUnsignedIntVector)
const
NekOptimize::GlobalOptParamSharedPtr
GetGlobalOptParam (void)
map< int, RobinBCInfoSharedPtrGetRobinBCInfo ()
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.
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.
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.
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.
void ExtractDataToCoeffs (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
 Extract the data in fielddata into the coeffs.
boost::shared_ptr< ExpListGetSharedThisPtr ()
 Returns a shared pointer to the current object.
boost::shared_ptr
< LibUtilities::SessionReader
GetSession ()
 Returns the session object.
boost::shared_ptr
< LibUtilities::Comm
GetComm ()
 Returns the comm object.
SpatialDomains::MeshGraphSharedPtr GetGraph ()
LibUtilities::BasisSharedPtr GetHomogeneousBasis (void)
boost::shared_ptr< ExpList > & GetPlane (int n)

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.
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.
void v_GetNormals (Array< OneD, Array< OneD, NekDouble > > &normals)
 Populate normals with the normals of all expansions.
- 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.
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.
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.
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.
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 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_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_SetUpPhysNormals ()
virtual void v_GetBoundaryToElmtMap (Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
virtual void v_ReadGlobalOptimizationParameters ()
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.
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::ofstream &outfile, std::string var="")
virtual void v_WriteTecplotZone (std::ofstream &outfile, int expansion)
virtual void v_WriteTecplotField (std::ofstream &outfile, int expansion)
virtual void v_WriteTecplotConnectivity (std::ofstream &outfile, int expansion)
virtual void v_WriteVtkPieceHeader (std::ofstream &outfile, int expansion)
virtual void v_WriteVtkPieceData (std::ofstream &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)
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.
virtual void v_ReadGlobalOptimizationParameters ()
virtual void v_SetUpPhysNormals ()
 Set up the normals on each expansion.
virtual void v_WriteVtkPieceHeader (std::ofstream &outfile, int expansion)
 const StdRegions::StdExpansionVector &locexp);

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.
LibUtilities::SessionReaderSharedPtr m_session
 Session.
SpatialDomains::MeshGraphSharedPtr m_graph
 Mesh associated with this expansion list.
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$.
int m_npoints
Array< OneD, NekDoublem_coeffs
 Concatenation of all local expansion coefficients.
Array< OneD, NekDoublem_phys
 The global expansion evaluated at the quadrature points.
bool m_physState
 The state of the array m_phys.
boost::shared_ptr
< LocalRegions::ExpansionVector
m_exp
 The list of local expansions.
Array< OneD, int > m_coeff_offset
 Offset of elemental data into the array m_coeffs.
Array< OneD, int > m_phys_offset
 Offset of elemental data into the array m_phys.
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]];.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
BlockMatrixMapShPtr m_blockMat
bool m_WaveSpace

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

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

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

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

:
ExpList(In,DeclareCoeffPhysArrays)
{
}
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 107 of file ExpList1D.cpp.

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

:
ExpList(pSession,graph1D)
{
int id=0;
const SpatialDomains::ExpansionMap &expansions
= graph1D->GetExpansions();
// For each element in the mesh, create a segment expansion using
// the supplied BasisKey and segment geometry.
SpatialDomains::ExpansionMap::const_iterator expIt;
for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
{
if ((SegmentGeom = boost::dynamic_pointer_cast<
expIt->second->m_geomShPtr)))
{
::AllocateSharedPtr(Ba,SegmentGeom);
seg->SetElmtId(id++);
(*m_exp).push_back(seg);
}
else
{
ASSERTL0(false,"dynamic cast to a SegGeom failed");
}
}
// Setup Default optimisation information.
int nel = GetExpSize();
// Allocate storage for data and populate element offset lists.
m_coeffs = Array<OneD, NekDouble>(m_ncoeffs);
m_phys = Array<OneD, NekDouble>(m_npoints);
}
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 176 of file ExpList1D.cpp.

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

:
ExpList(pSession,graph1D)
{
int id=0;
// Retrieve the list of expansions
const SpatialDomains::ExpansionMap &expansions
= graph1D->GetExpansions();
// Process each expansion in the graph
SpatialDomains::ExpansionMap::const_iterator expIt;
for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
{
// Retrieve basis key from expansion
LibUtilities::BasisKey bkey = expIt->second->m_basisKeyVector[0];
if ((SegmentGeom = boost::dynamic_pointer_cast<
expIt->second->m_geomShPtr)))
{
::AllocateSharedPtr(bkey, SegmentGeom);
// Assign next ID
seg->SetElmtId(id++);
// Add the expansion
(*m_exp).push_back(seg);
}
else
{
ASSERTL0(false,"dynamic cast to a SegGeom failed");
}
}
// Setup Default optimisation information.
int nel = GetExpSize();
// set up offset arrays.
if(DeclareCoeffPhysArrays)
{
// Set up m_coeffs, m_phys.
m_coeffs = Array<OneD, NekDouble>(m_ncoeffs);
m_phys = Array<OneD, NekDouble>(m_npoints);
}
}
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 258 of file ExpList1D.cpp.

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

:
ExpList(pSession,graph1D)
{
int j,id=0;
SpatialDomains::CompositeMap::const_iterator compIt;
// Retrieve the list of expansions
const SpatialDomains::ExpansionMap &expansions
= graph1D->GetExpansions(var);
// Process each composite region in domain
for(compIt = domain.begin(); compIt != domain.end(); ++compIt)
{
comp = compIt->second;
// Process each expansion in the graph
for(j = 0; j < compIt->second->size(); ++j)
{
SpatialDomains::ExpansionMap::const_iterator expIt;
if((SegmentGeom = boost::dynamic_pointer_cast<
(*compIt->second)[j])))
{
// Retrieve basis key from expansion and define expansion
if((expIt = expansions.find(SegmentGeom->GetGlobalID())) != expansions.end())
{
LibUtilities::BasisKey bkey = expIt->second->m_basisKeyVector[0];
if(SetToOneSpaceDimension)
{
SegmentGeom->GenerateOneSpaceDimGeom();
::AllocateSharedPtr(bkey, OneDSegmentGeom);
}
else
{
::AllocateSharedPtr(bkey, SegmentGeom);
}
}
else
{
ASSERTL0(false,"Failed to find basis key");
}
}
else
{
ASSERTL0(false,"Failed to dynamic cast geometry to SegGeom");
}
// Assign next ID
seg->SetElmtId(id++);
// Add the expansion
(*m_exp).push_back(seg);
}
}
// Setup Default optimisation information.
int nel = GetExpSize();
// set up offset arrays.
if(DeclareCoeffPhysArrays)
{
// Set up m_coeffs, m_phys.
m_coeffs = Array<OneD, NekDouble>(m_ncoeffs);
m_phys = Array<OneD, NekDouble>(m_npoints);
}
}
Nektar::MultiRegions::ExpList1D::ExpList1D ( 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 361 of file ExpList1D.cpp.

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

:
{
int j, id=0;
SpatialDomains::CompositeMap::const_iterator compIt;
// Process each composite region.
for(compIt = domain.begin(); compIt != domain.end(); ++compIt)
{
comp = compIt->second;
// Process each expansion in the region.
for(j = 0; j < compIt->second->size(); ++j)
{
if((SegmentGeom = boost::dynamic_pointer_cast<
(*compIt->second)[j])))
{
// Retrieve the basis key from the expansion.
= boost::dynamic_pointer_cast<SpatialDomains::MeshGraph2D>(graph2D)->GetEdgeBasisKey(SegmentGeom, variable);
::AllocateSharedPtr(bkey, SegmentGeom);
// Add the segment to the expansion list.
seg->SetElmtId(id++);
(*m_exp).push_back(seg);
}
else
{
ASSERTL0(false,"dynamic cast to a SegGeom failed");
}
}
}
// Setup Default optimisation information.
int nel = GetExpSize();
// Allocate storage for data and populate element offset lists.
// Set up m_coeffs, m_phys.
if(DeclareCoeffPhysArrays)
{
m_coeffs = Array<OneD, NekDouble>(m_ncoeffs);
m_phys = Array<OneD, NekDouble>(m_npoints);
}
}
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 437 of file ExpList1D.cpp.

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

:
{
int i, j, id, elmtid = 0;
set<int> edgesDone;
map<int,int> EdgeDone;
map<int,int> NormalSet;
// First loop over boundary conditions to renumber
// Dirichlet boundaries
for(i = 0; i < bndCond.num_elements(); ++i)
{
if(bndCond[i]->GetBoundaryConditionType()
{
for(j = 0; j < bndConstraint[i]->GetExpSize(); ++j)
{
LibUtilities::BasisKey bkey = bndConstraint[i]
->GetExp(j)->GetBasis(0)->GetBasisKey();
exp1D = bndConstraint[i]->GetExp(j)->
as<LocalRegions::Expansion1D>();
segGeom = exp1D->GetGeom1D();
::AllocateSharedPtr(bkey, segGeom);
edgesDone.insert(segGeom->GetEid());
seg->SetElmtId(elmtid++);
(*m_exp).push_back(seg);
}
}
}
LibUtilities::BasisKey> > edgeOrders;
for(i = 0; i < locexp.size(); ++i)
{
exp2D = locexp[i]->as<LocalRegions::Expansion2D>();
for(j = 0; j < locexp[i]->GetNedges(); ++j)
{
segGeom = exp2D->GetGeom2D()->GetEdge(j);
id = segGeom->GetEid();
// Ignore Dirichlet edges
if (edgesDone.count(id) != 0)
{
continue;
}
it = edgeOrders.find(id);
if (it == edgeOrders.end())
{
edgeOrders.insert(std::make_pair(id, std::make_pair(
segGeom, locexp[i]->DetEdgeBasisKey(j))));
}
else // variable modes/points
{
= locexp[i]->DetEdgeBasisKey(j);
= it->second.second;
int np1 = edge .GetNumPoints();
int np2 = existing.GetNumPoints();
int nm1 = edge .GetNumModes ();
int nm2 = existing.GetNumModes ();
if (np2 >= np1 && nm2 >= nm1)
{
continue;
}
else if (np2 < np1 && nm2 < nm1)
{
it->second.second = edge;
}
else
{
ASSERTL0(false,
"inappropriate number of points/modes (max "
"num of points is not set with max order)");
}
}
}
}
pSession->GetComm()->GetRowComm();
int nproc = vComm->GetSize(); // number of processors
int edgepr = vComm->GetRank(); // ID processor
if (nproc > 1)
{
int eCnt = 0;
// Count the number of edges on each partition
for(i = 0; i < locexp.size(); ++i)
{
eCnt += locexp[i]->GetNedges();
}
// Set up the offset and the array that will contain the list of
// edge IDs, then reduce this across processors.
Array<OneD, int> edgesCnt(nproc, 0);
edgesCnt[edgepr] = eCnt;
vComm->AllReduce(edgesCnt, LibUtilities::ReduceSum);
// Set up offset array.
int totEdgeCnt = Vmath::Vsum(nproc, edgesCnt, 1);
Array<OneD, int> eTotOffsets(nproc,0);
for (i = 1; i < nproc; ++i)
{
eTotOffsets[i] = eTotOffsets[i-1] + edgesCnt[i-1];
}
// Local list of the edges per element
Array<OneD, int> EdgesTotID(totEdgeCnt, 0);
Array<OneD, int> EdgesTotNm(totEdgeCnt, 0);
Array<OneD, int> EdgesTotPnts(totEdgeCnt, 0);
int cntr = eTotOffsets[edgepr];
for(i = 0; i < locexp.size(); ++i)
{
exp2D = locexp[i]->as<LocalRegions::Expansion2D>();
int nedges = locexp[i]->GetNedges();
for(j = 0; j < nedges; ++j, ++cntr)
{
locexp[i]->DetEdgeBasisKey(j);
EdgesTotID [cntr] = exp2D->GetGeom2D()->GetEid(j);
EdgesTotNm [cntr] = bkeyEdge.GetNumModes();
EdgesTotPnts[cntr] = bkeyEdge.GetNumPoints();
}
}
vComm->AllReduce(EdgesTotID, LibUtilities::ReduceSum);
vComm->AllReduce(EdgesTotNm, LibUtilities::ReduceSum);
vComm->AllReduce(EdgesTotPnts, LibUtilities::ReduceSum);
for (i = 0; i < totEdgeCnt; ++i)
{
it = edgeOrders.find(EdgesTotID[i]);
if (it == edgeOrders.end())
{
continue;
}
= it->second.second;
existing.GetBasisType(), EdgesTotNm[i],
LibUtilities::PointsKey(EdgesTotPnts[i],
existing.GetPointsType()));
int np1 = edge .GetNumPoints();
int np2 = existing.GetNumPoints();
int nm1 = edge .GetNumModes ();
int nm2 = existing.GetNumModes ();
if (np2 >= np1 && nm2 >= nm1)
{
continue;
}
else if (np2 < np1 && nm2 < nm1)
{
it->second.second = edge;
}
else
{
ASSERTL0(false,
"inappropriate number of points/modes (max "
"num of points is not set with max order)");
}
}
}
for (it = edgeOrders.begin(); it != edgeOrders.end(); ++it)
{
::AllocateSharedPtr(it->second.second, it->second.first);
seg->SetElmtId(elmtid++);
(*m_exp).push_back(seg);
}
// Setup Default optimisation information.
int nel = GetExpSize();
// Set up offset information and array sizes
// Set up m_coeffs, m_phys.
if(DeclareCoeffPhysArrays)
{
m_coeffs = Array<OneD, NekDouble>(m_ncoeffs);
m_phys = Array<OneD, NekDouble>(m_npoints);
}
}
Nektar::MultiRegions::ExpList1D::~ExpList1D ( )
virtual

Destructor.

Definition at line 699 of file ExpList1D.cpp.

{
}

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

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

Referenced by PostProcess().

{
int i,j,r;
// Get the number of elements in the domain
int num_elm = GetExpSize();
// initializing the outarray
for(i = 0; i < outarray.num_elements(); i++)
{
outarray[i] = 0.0;
}
// Make a copy for further modification
int x_size = inarray2.num_elements();
Array<OneD,NekDouble> x_values_cp(x_size);
// Determining the element to which the x belongs
Array<OneD,int> x_elm(x_size);
for(i = 0; i < x_size; i++ )
{
x_elm[i] = (int)floor(inarray2[i]/h);
}
// Clamp indices periodically
for(i = 0; i < x_size; i++)
{
while(x_elm[i] < 0)
{
x_elm[i] += num_elm;
}
while(x_elm[i] >= num_elm)
{
x_elm[i] -= num_elm ;
}
}
// Map the values of x to [-1 1] on its interval
for(i = 0; i < x_size; i++)
{
x_values_cp[i] = (inarray2[i]/h - floor(inarray2[i]/h))*2 - 1.0;
}
// Evaluate the jocobi polynomials
// (Evaluating the base at some points other than the quadrature
// points). Should it be added to the base class????
Array<TwoD,NekDouble> jacobi_poly(nmodes,x_size);
for(i = 0; i < nmodes; i++)
{
Polylib::jacobfd(x_size,x_values_cp.get(),
jacobi_poly.get()+i*x_size,NULL,i,0.0,0.0);
}
// Evaluate the function values
for(r = 0; r < nmodes; r++)
{
for(j = 0; j < x_size; j++)
{
int index = ((x_elm[j])*nmodes)+r;
outarray[j] += inarray1[index]*jacobi_poly[r][j];
}
}
}
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 713 of file ExpList1D.cpp.

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

{
int i,j,r;
// get the local element expansion of the elmId element
// Get the quadrature points and weights required for integration
int quad_npoints = elmExp->GetTotPoints();
LibUtilities::PointsKey quadPointsKey(quad_npoints,
elmExp->GetPointsType(0));
Array<OneD,NekDouble> quad_points
= LibUtilities::PointsManager()[quadPointsKey]->GetZ();
Array<OneD,NekDouble> quad_weights
= LibUtilities::PointsManager()[quadPointsKey]->GetW();
// Declare variable for the local kernel breaks
int kernel_width = kernel->GetKernelWidth();
Array<OneD,NekDouble> local_kernel_breaks(kernel_width+1);
// Declare variable for the transformed quadrature points
Array<OneD,NekDouble> mapped_quad_points(quad_npoints);
// For each evaluation point
for(i = 0; i < inarray.num_elements(); i++)
{
// Move the center of the kernel to the current point
kernel->MoveKernelCenter(inarray[i],local_kernel_breaks);
// Find the mesh breaks under the kernel support
Array<OneD,NekDouble> mesh_breaks;
kernel->FindMeshUnderKernel(local_kernel_breaks,h,mesh_breaks);
// Sort the total breaks for integration purposes
int total_nbreaks = local_kernel_breaks.num_elements() +
mesh_breaks.num_elements();
// number of the total breaks
Array<OneD,NekDouble> total_breaks(total_nbreaks);
kernel->Sort(local_kernel_breaks,mesh_breaks,total_breaks);
// Integrate the product of kernel and function over the total
// breaks
NekDouble integral_value = 0.0;
for(j = 0; j < total_breaks.num_elements()-1; j++)
{
NekDouble a = total_breaks[j];
NekDouble b = total_breaks[j+1];
// Map the quadrature points to the appropriate interval
for(r = 0; r < quad_points.num_elements(); r++)
{
mapped_quad_points[r]
= (quad_points[r] + 1.0) * 0.5 * (b - a) + a;
}
// Evaluate the function at the transformed quadrature
// points
Array<OneD,NekDouble> u_value(quad_npoints);
Array<OneD,NekDouble> coeffs = GetCoeffs();
PeriodicEval(coeffs,mapped_quad_points,h,
elmExp->GetBasisNumModes(0),u_value);
// Evaluate the kernel at the transformed quadrature points
Array<OneD,NekDouble> k_value(quad_npoints);
kernel->EvaluateKernel(mapped_quad_points,h,k_value);
// Integrate
for(r = 0; r < quad_npoints; r++)
{
integral_value += (b - a) * 0.5 * k_value[r]
* u_value[r] * quad_weights[r];
}
}
outarray[i] = integral_value/h;
}
}
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 675 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().

{
int i;
// Set up offset information and array sizes
m_coeff_offset = Array<OneD,int>(m_exp->size());
m_phys_offset = Array<OneD,int>(m_exp->size());
m_offset_elmt_id = Array<OneD,int>(m_exp->size());
for(i = 0; i < m_exp->size(); ++i)
{
m_phys_offset [i] = m_npoints;
m_offset_elmt_id[i] = i;
m_ncoeffs += (*m_exp)[i]->GetNcoeffs();
m_npoints += (*m_exp)[i]->GetTotPoints();
}
}
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 1040 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.

{
int i,j,k,e_npoints,offset;
Array<OneD,Array<OneD,NekDouble> > locnormals;
Array<OneD,Array<OneD,NekDouble> > locnormals2;
Array<OneD,Array<OneD,NekDouble> > Norms;
// Assume whole array is of same coordinate dimension
int coordim = GetCoordim(0);
ASSERTL1(normals.num_elements() >= coordim,
"Output vector does not have sufficient dimensions to "
"match coordim");
for (i = 0; i < m_exp->size(); ++i)
{
int edgeNumber = loc_exp->GetLeftAdjacentElementEdge();
// Get the number of points and normals for this expansion.
e_npoints = (*m_exp)[i]->GetNumPoints(0);
locnormals = loc_elmt->GetEdgeNormal(edgeNumber);
int e_nmodes = loc_exp->GetBasis(0)->GetNumModes();
int loc_nmodes = loc_elmt->GetBasis(0)->GetNumModes();
if (e_nmodes != loc_nmodes)
{
if (loc_exp->GetRightAdjacentElementEdge() >= 0)
{
loc_exp->GetRightAdjacentElementExp();
int EdgeNumber = loc_exp->GetRightAdjacentElementEdge();
// Serial case: right element is connected so we can
// just grab that normal.
locnormals = loc_elmt->GetEdgeNormal(EdgeNumber);
offset = m_phys_offset[i];
// Process each point in the expansion.
for (j = 0; j < e_npoints; ++j)
{
// Process each spatial dimension and copy the values
// into the output array.
for (k = 0; k < coordim; ++k)
{
normals[k][offset + j] = -locnormals[k][j];
}
}
}
else
{
// Parallel case: need to interpolate normal.
Array<OneD, Array<OneD, NekDouble> > normal(coordim);
for (int p = 0; p < coordim; ++p)
{
normal[p] = Array<OneD, NekDouble>(e_npoints,0.0);
loc_exp->GetBasis(0)->GetPointsKey();
loc_elmt->GetBasis(0)->GetPointsKey();
locnormals[p],
to_key,
normal[p]);
}
offset = m_phys_offset[i];
// Process each point in the expansion.
for (j = 0; j < e_npoints; ++j)
{
// Process each spatial dimension and copy the values
// into the output array.
for (k = 0; k < coordim; ++k)
{
normals[k][offset + j] = normal[k][j];
}
}
}
}
else
{
// Get the physical data offset for this expansion.
offset = m_phys_offset[i];
// Process each point in the expansion.
for (j = 0; j < e_npoints; ++j)
{
// Process each spatial dimension and copy the values
// into the output array.
for (k = 0; k < coordim; ++k)
{
normals[k][offset + j] = locnormals[k][j];
}
}
}
}
}
void Nektar::MultiRegions::ExpList1D::v_ReadGlobalOptimizationParameters ( )
privatevirtual

Definition at line 1150 of file ExpList1D.cpp.

{
// Array<OneD, int> NumShape(1,0);
// NumShape[0] = GetExpSize();
//
// int one = 1;
// m_globalOptParam = MemoryManager<NekOptimize::GlobalOptParam>
// ::AllocateSharedPtr(m_session,one,NumShape);
}
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.

Definition at line 915 of file ExpList1D.cpp.

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

{
int i, j;
for (i = 0; i < m_exp->size(); ++i)
{
for (j = 0; j < (*m_exp)[i]->GetNverts(); ++j)
{
(*m_exp)[i]->ComputeVertexNormal(j);
}
}
}
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 937 of file ExpList1D.cpp.

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

{
int i,j,k,e_npoints,offset;
Array<OneD,NekDouble> normals;
// Assume whole array is of same coordimate dimension
int coordim = GetCoordim(0);
ASSERTL1(Vec.num_elements() >= coordim,
"Input vector does not have sufficient dimensions to "
"match coordim");
// Process each expansion
for(i = 0; i < m_exp->size(); ++i)
{
// Get the number of points in the expansion and the normals.
e_npoints = (*m_exp)[i]->GetNumPoints(0);
normals = (*m_exp)[i]->GetPhysNormals();
// Get the physical data offset of the expansion in m_phys.
offset = m_phys_offset[i];
// Compute each data point.
for(j = 0; j < e_npoints; ++j)
{
// Calculate normal velocity.
Vn = 0.0;
for(k = 0; k < coordim; ++k)
{
Vn += Vec[k][offset+j]*normals[k*e_npoints + j];
}
// Upwind based on direction of normal velocity.
if(Vn > 0.0)
{
Upwind[offset + j] = Fwd[offset + j];
}
else
{
Upwind[offset + j] = Bwd[offset + j];
}
}
}
}
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 1000 of file ExpList1D.cpp.

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

{
int i,j,e_npoints,offset;
Array<OneD,NekDouble> normals;
// Process each expansion.
for(i = 0; i < m_exp->size(); ++i)
{
// Get the number of points and the data offset.
e_npoints = (*m_exp)[i]->GetNumPoints(0);
offset = m_phys_offset[i];
// Process each point in the expansion.
for(j = 0; j < e_npoints; ++j)
{
// Upwind based on one-dimensional velocity.
if(Vn[offset + j] > 0.0)
{
Upwind[offset + j] = Fwd[offset + j];
}
else
{
Upwind[offset + j] = Bwd[offset + j];
}
}
}
}
void Nektar::MultiRegions::ExpList1D::v_WriteVtkPieceHeader ( std::ofstream &  outfile,
int  expansion 
)
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.

Definition at line 1169 of file ExpList1D.cpp.

{
int i,j;
int nquad0 = (*m_exp)[expansion]->GetNumPoints(0);
int ntot = nquad0;
int ntotminus = (nquad0-1);
Array<OneD,NekDouble> coords[3];
coords[0] = Array<OneD,NekDouble>(ntot);
coords[1] = Array<OneD,NekDouble>(ntot);
coords[2] = Array<OneD,NekDouble>(ntot);
(*m_exp)[expansion]->GetCoords(coords[0],coords[1],coords[2]);
outfile << " <Piece NumberOfPoints=\""
<< ntot << "\" NumberOfCells=\""
<< ntotminus << "\">" << endl;
outfile << " <Points>" << endl;
outfile << " <DataArray type=\"Float32\" "
<< "NumberOfComponents=\"3\" format=\"ascii\">" << endl;
outfile << " ";
for (i = 0; i < ntot; ++i)
{
for (j = 0; j < 3; ++j)
{
outfile << setprecision(8) << scientific
<< (float)coords[j][i] << " ";
}
outfile << endl;
}
outfile << endl;
outfile << " </DataArray>" << endl;
outfile << " </Points>" << endl;
outfile << " <Cells>" << endl;
outfile << " <DataArray type=\"Int32\" "
<< "Name=\"connectivity\" format=\"ascii\">" << endl;
for (i = 0; i < nquad0-1; ++i)
{
outfile << i << " " << i+1 << endl;
}
outfile << endl;
outfile << " </DataArray>" << endl;
outfile << " <DataArray type=\"Int32\" "
<< "Name=\"offsets\" format=\"ascii\">" << endl;
for (i = 0; i < ntotminus; ++i)
{
outfile << i*2+2 << " ";
}
outfile << endl;
outfile << " </DataArray>" << endl;
outfile << " <DataArray type=\"UInt8\" "
<< "Name=\"types\" format=\"ascii\">" << endl;
for (i = 0; i < ntotminus; ++i)
{
outfile << "3 ";
}
outfile << endl;
outfile << " </DataArray>" << endl;
outfile << " </Cells>" << endl;
outfile << " <PointData>" << endl;
}

Member Data Documentation

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

Definition at line 184 of file ExpList1D.h.

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

Definition at line 186 of file ExpList1D.h.