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

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

75  :
76  ExpList()
77  {
78  SetExpType(e1D);
79  }
ExpList()
The default constructor.
Definition: ExpList.cpp:93
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:251
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().

85  :
86  ExpList(In,DeclareCoeffPhysArrays)
87  {
88  SetExpType(e1D);
89  }
ExpList()
The default constructor.
Definition: ExpList.cpp:93
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:251
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 94 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().

96  :
97  ExpList(In, eIDs, DeclareCoeffPhysArrays)
98  {
99  SetExpType(e1D);
100 
101  // Setup Default optimisation information.
102  int nel = GetExpSize();
105 
106  // Allocate storage for data and populate element offset lists.
108 
111  }
ExpList()
The default constructor.
Definition: ExpList.cpp:93
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:1001
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1896
void SetCoeffPhysOffsets(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
Definition: ExpList1D.cpp:707
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:2950
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:251
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 129 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().

131  :
132  ExpList(pSession,graph1D)
133  {
134  SetExpType(e1D);
135 
136  int id=0;
139 
140  const SpatialDomains::ExpansionMap &expansions
141  = graph1D->GetExpansions();
142 
143  // For each element in the mesh, create a segment expansion using
144  // the supplied BasisKey and segment geometry.
145  SpatialDomains::ExpansionMap::const_iterator expIt;
146  for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
147  {
148  if ((SegmentGeom = boost::dynamic_pointer_cast<
150  expIt->second->m_geomShPtr)))
151  {
153  ::AllocateSharedPtr(Ba,SegmentGeom);
154  seg->SetElmtId(id++);
155  (*m_exp).push_back(seg);
156  }
157  else
158  {
159  ASSERTL0(false,"dynamic cast to a SegGeom failed");
160  }
161  }
162 
163  // Setup Default optimisation information.
164  int nel = GetExpSize();
167 
168  // Allocate storage for data and populate element offset lists.
170 
173 
176  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
ExpList()
The default constructor.
Definition: ExpList.cpp:93
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h: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:1896
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:707
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:2950
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:251
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 199 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().

201  :
202  ExpList(pSession,graph1D)
203  {
204  SetExpType(e1D);
205 
206  int id=0;
209 
210  // Retrieve the list of expansions
211  const SpatialDomains::ExpansionMap &expansions
212  = graph1D->GetExpansions();
213 
214  // Process each expansion in the graph
215  SpatialDomains::ExpansionMap::const_iterator expIt;
216  for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
217  {
218  // Retrieve basis key from expansion
219  LibUtilities::BasisKey bkey = expIt->second->m_basisKeyVector[0];
220 
221  if ((SegmentGeom = boost::dynamic_pointer_cast<
223  expIt->second->m_geomShPtr)))
224  {
226  ::AllocateSharedPtr(bkey, SegmentGeom);
227 
228  // Assign next ID
229  seg->SetElmtId(id++);
230 
231  // Add the expansion
232  (*m_exp).push_back(seg);
233  }
234  else
235  {
236  ASSERTL0(false,"dynamic cast to a SegGeom failed");
237  }
238  }
239 
240  // Setup Default optimisation information.
241  int nel = GetExpSize();
244 
245  // set up offset arrays.
247 
248  if(DeclareCoeffPhysArrays)
249  {
250  // Set up m_coeffs, m_phys.
253  }
254 
257  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
ExpList()
The default constructor.
Definition: ExpList.cpp:93
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h: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:1896
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:707
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:2950
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:251
Describes the specification for a Basis.
Definition: Basis.h:50
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 282 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().

287  :
288  ExpList(pSession,graph1D)
289  {
290  int j,id=0;
294  SpatialDomains::CompositeMap::const_iterator compIt;
295 
296  // Retrieve the list of expansions
297  const SpatialDomains::ExpansionMap &expansions
298  = graph1D->GetExpansions(var);
299 
300  // Process each composite region in domain
301  for(compIt = domain.begin(); compIt != domain.end(); ++compIt)
302  {
303  comp = compIt->second;
304 
305  // Process each expansion in the graph
306  for(j = 0; j < compIt->second->size(); ++j)
307  {
308  SpatialDomains::ExpansionMap::const_iterator expIt;
309 
310  if((SegmentGeom = boost::dynamic_pointer_cast<
312  (*compIt->second)[j])))
313  {
314  // Retrieve basis key from expansion and define expansion
315  if((expIt = expansions.find(SegmentGeom->GetGlobalID())) != expansions.end())
316  {
317  LibUtilities::BasisKey bkey = expIt->second->m_basisKeyVector[0];
318 
319  if(SetToOneSpaceDimension)
320  {
321  SpatialDomains::SegGeomSharedPtr OneDSegmentGeom =
322  SegmentGeom->GenerateOneSpaceDimGeom();
323 
325  ::AllocateSharedPtr(bkey, OneDSegmentGeom);
326  }
327  else
328  {
330  ::AllocateSharedPtr(bkey, SegmentGeom);
331  }
332  }
333  else
334  {
335  ASSERTL0(false,"Failed to find basis key");
336  }
337 
338  }
339  else
340  {
341  ASSERTL0(false,"Failed to dynamic cast geometry to SegGeom");
342  }
343 
344 
345  // Assign next ID
346  seg->SetElmtId(id++);
347 
348  // Add the expansion
349  (*m_exp).push_back(seg);
350  }
351  }
352 
353  // Setup Default optimisation information.
354  int nel = GetExpSize();
357 
358  // set up offset arrays.
360 
361  if(DeclareCoeffPhysArrays)
362  {
363  // Set up m_coeffs, m_phys.
366  }
367 
370  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
ExpList()
The default constructor.
Definition: ExpList.cpp:93
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h: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:1896
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:707
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:2950
Describes the specification for a Basis.
Definition: Basis.h:50
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 386 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().

390  :
391  ExpList(pSession,graph2D)
392  {
393  SetExpType(e1D);
394 
395  m_graph = graph2D;
396 
397  int j, id=0;
399  SpatialDomains::CompositeMap::const_iterator compIt;
402 
403  // Process each composite region.
404  for(compIt = domain.begin(); compIt != domain.end(); ++compIt)
405  {
406  comp = compIt->second;
407  // Process each expansion in the region.
408  for(j = 0; j < compIt->second->size(); ++j)
409  {
410  if((SegmentGeom = boost::dynamic_pointer_cast<
412  (*compIt->second)[j])))
413  {
414  // Retrieve the basis key from the expansion.
416  = boost::dynamic_pointer_cast<SpatialDomains::MeshGraph2D>(graph2D)->GetEdgeBasisKey(SegmentGeom, variable);
417 
419  ::AllocateSharedPtr(bkey, SegmentGeom);
420 
421  // Add the segment to the expansion list.
422  seg->SetElmtId(id++);
423  (*m_exp).push_back(seg);
424  }
425  else
426  {
427  ASSERTL0(false,"dynamic cast to a SegGeom failed");
428  }
429  }
430 
431  }
432 
433  // Setup Default optimisation information.
434  int nel = GetExpSize();
437 
438  // Allocate storage for data and populate element offset lists.
440 
441  // Set up m_coeffs, m_phys.
442  if(DeclareCoeffPhysArrays)
443  {
446  }
447 
449  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
ExpList()
The default constructor.
Definition: ExpList.cpp:93
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h: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:1896
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:707
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:2950
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:251
Describes the specification for a Basis.
Definition: Basis.h:50
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 467 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().

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

Destructor.

Definition at line 731 of file ExpList1D.cpp.

732  {
733  }

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

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

Referenced by PostProcess().

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

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

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

708  {
709  int i;
710 
711  // Set up offset information and array sizes
715 
716  m_ncoeffs = m_npoints = 0;
717 
718  for(i = 0; i < m_exp->size(); ++i)
719  {
721  m_phys_offset [i] = m_npoints;
722  m_offset_elmt_id[i] = i;
723  m_ncoeffs += (*m_exp)[i]->GetNcoeffs();
724  m_npoints += (*m_exp)[i]->GetTotPoints();
725  }
726  }
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 1072 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.

1074  {
1075  int i,j,k,e_npoints,offset;
1077  Array<OneD,Array<OneD,NekDouble> > locnormals;
1078  Array<OneD,Array<OneD,NekDouble> > locnormals2;
1080  // Assume whole array is of same coordinate dimension
1081  int coordim = GetCoordim(0);
1082 
1083  ASSERTL1(normals.num_elements() >= coordim,
1084  "Output vector does not have sufficient dimensions to "
1085  "match coordim");
1086 
1087  for (i = 0; i < m_exp->size(); ++i)
1088  {
1090 
1092  loc_exp->GetLeftAdjacentElementExp();
1093 
1094  int edgeNumber = loc_exp->GetLeftAdjacentElementEdge();
1095 
1096  // Get the number of points and normals for this expansion.
1097  e_npoints = (*m_exp)[i]->GetNumPoints(0);
1098 
1099  locnormals = loc_elmt->GetEdgeNormal(edgeNumber);
1100  int e_nmodes = loc_exp->GetBasis(0)->GetNumModes();
1101  int loc_nmodes = loc_elmt->GetBasis(0)->GetNumModes();
1102 
1103  if (e_nmodes != loc_nmodes)
1104  {
1105  if (loc_exp->GetRightAdjacentElementEdge() >= 0)
1106  {
1108  loc_exp->GetRightAdjacentElementExp();
1109 
1110  int EdgeNumber = loc_exp->GetRightAdjacentElementEdge();
1111  // Serial case: right element is connected so we can
1112  // just grab that normal.
1113  locnormals = loc_elmt->GetEdgeNormal(EdgeNumber);
1114 
1115  offset = m_phys_offset[i];
1116 
1117  // Process each point in the expansion.
1118  for (j = 0; j < e_npoints; ++j)
1119  {
1120  // Process each spatial dimension and copy the values
1121  // into the output array.
1122  for (k = 0; k < coordim; ++k)
1123  {
1124  normals[k][offset + j] = -locnormals[k][j];
1125  }
1126  }
1127  }
1128  else
1129  {
1130  // Parallel case: need to interpolate normal.
1131  Array<OneD, Array<OneD, NekDouble> > normal(coordim);
1132 
1133  for (int p = 0; p < coordim; ++p)
1134  {
1135  normal[p] = Array<OneD, NekDouble>(e_npoints,0.0);
1136  LibUtilities::PointsKey to_key =
1137  loc_exp->GetBasis(0)->GetPointsKey();
1138  LibUtilities::PointsKey from_key =
1139  loc_elmt->GetBasis(0)->GetPointsKey();
1140  LibUtilities::Interp1D(from_key,
1141  locnormals[p],
1142  to_key,
1143  normal[p]);
1144  }
1145 
1146  offset = m_phys_offset[i];
1147 
1148  // Process each point in the expansion.
1149  for (j = 0; j < e_npoints; ++j)
1150  {
1151  // Process each spatial dimension and copy the values
1152  // into the output array.
1153  for (k = 0; k < coordim; ++k)
1154  {
1155  normals[k][offset + j] = normal[k][j];
1156  }
1157  }
1158  }
1159  }
1160  else
1161  {
1162  // Get the physical data offset for this expansion.
1163  offset = m_phys_offset[i];
1164 
1165  // Process each point in the expansion.
1166  for (j = 0; j < e_npoints; ++j)
1167  {
1168  // Process each spatial dimension and copy the values
1169  // into the output array.
1170  for (k = 0; k < coordim; ++k)
1171  {
1172  normals[k][offset + j] = locnormals[k][j];
1173  }
1174  }
1175  }
1176  }
1177  }
Expansion2DSharedPtr GetLeftAdjacentElementExp() const
Definition: Expansion1D.h:123
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
Defines a specification for a set of points.
Definition: Points.h:58
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:1794
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
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 1182 of file ExpList1D.cpp.

1183  {
1184 // Array<OneD, int> NumShape(1,0);
1185 // NumShape[0] = GetExpSize();
1186 //
1187 // int one = 1;
1188 // m_globalOptParam = MemoryManager<NekOptimize::GlobalOptParam>
1189 // ::AllocateSharedPtr(m_session,one,NumShape);
1190  }
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 947 of file ExpList1D.cpp.

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

948  {
949  int i, j;
950  for (i = 0; i < m_exp->size(); ++i)
951  {
952  for (j = 0; j < (*m_exp)[i]->GetNverts(); ++j)
953  {
954  (*m_exp)[i]->ComputeVertexNormal(j);
955  }
956  }
957  }
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 969 of file ExpList1D.cpp.

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

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

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

1037  {
1038  int i,j,e_npoints,offset;
1039  Array<OneD,NekDouble> normals;
1040 
1041  // Process each expansion.
1042  for(i = 0; i < m_exp->size(); ++i)
1043  {
1044  // Get the number of points and the data offset.
1045  e_npoints = (*m_exp)[i]->GetNumPoints(0);
1046  offset = m_phys_offset[i];
1047 
1048  // Process each point in the expansion.
1049  for(j = 0; j < e_npoints; ++j)
1050  {
1051  // Upwind based on one-dimensional velocity.
1052  if(Vn[offset + j] > 0.0)
1053  {
1054  Upwind[offset + j] = Fwd[offset + j];
1055  }
1056  else
1057  {
1058  Upwind[offset + j] = Bwd[offset + j];
1059  }
1060  }
1061  }
1062  }
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 1201 of file ExpList1D.cpp.

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

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.