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

Protected Member Functions

void v_Upwind (const Array< OneD, const 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 void v_ExtractTracePhys (Array< OneD, NekDouble > &outarray)
 
virtual void v_ExtractTracePhys (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_MultiplyByInvMassMatrix (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
 
virtual void v_HelmSolve (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const FlagList &flags, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff, const Array< OneD, const NekDouble > &dirForcing)
 
virtual void v_LinearAdvectionDiffusionReactionSolve (const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 
virtual void v_LinearAdvectionReactionSolve (const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 
virtual void v_ImposeDirichletConditions (Array< OneD, NekDouble > &outarray)
 
virtual void v_FillBndCondFromField ()
 
virtual void v_Reset ()
 Reset geometry information, metrics, matrix managers and geometry information. More...
 
virtual void v_LocalToGlobal (void)
 
virtual void v_GlobalToLocal (void)
 
virtual void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
 
virtual void v_BwdTrans_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
 
virtual void v_FwdTrans_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_SmoothField (Array< OneD, NekDouble > &field)
 
virtual void v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
 
virtual void v_IProductWRTBase_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_GeneralMatrixOp (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
 
virtual void v_GetCoords (Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
 
virtual void v_PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
 
virtual void v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
 
virtual void v_PhysDeriv (Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
 
virtual void v_HomogeneousFwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
 
virtual void v_HomogeneousBwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
 
virtual void v_DealiasedProd (const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
 
virtual void v_GetBCValues (Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
 
virtual void v_NormVectorIProductWRTBase (Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
 
virtual void v_GetBoundaryToElmtMap (Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
 
virtual std::vector< LibUtilities::FieldDefinitionsSharedPtrv_GetFieldDefinitions (void)
 
virtual void v_GetFieldDefinitions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef)
 
virtual void v_AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
 
virtual void v_AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, Array< OneD, NekDouble > &coeffs)
 
virtual void v_ExtractDataToCoeffs (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
 Extract data from raw field data into expansion list. More...
 
virtual void v_ExtractCoeffsToCoeffs (const boost::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
 
virtual void v_WriteTecplotHeader (std::ostream &outfile, std::string var="")
 
virtual void v_WriteTecplotZone (std::ostream &outfile, int expansion)
 
virtual void v_WriteTecplotField (std::ostream &outfile, int expansion)
 
virtual void v_WriteTecplotConnectivity (std::ostream &outfile, int expansion)
 
virtual void v_WriteVtkPieceData (std::ostream &outfile, int expansion, std::string var)
 
virtual NekDouble v_L2 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
 
virtual NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray)
 
virtual Array< OneD, const NekDoublev_HomogeneousEnergy (void)
 
virtual LibUtilities::TranspositionSharedPtr v_GetTransposition (void)
 
virtual NekDouble v_GetHomoLen (void)
 
virtual Array< OneD, const unsigned int > v_GetZIDs (void)
 
virtual Array< OneD, const unsigned int > v_GetYIDs (void)
 
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. 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::ExpansionVectorm_exp
 The list of local expansions. More...
 
Collections::CollectionVector m_collections
 
std::vector< int > m_coll_coeff_offset
 Offset of elemental data into the array m_coeffs. More...
 
std::vector< int > m_coll_phys_offset
 Offset of elemental data into the array m_phys. More...
 
Array< OneD, int > m_coeff_offset
 Offset of elemental data into the array m_coeffs. More...
 
Array< OneD, int > m_phys_offset
 Offset of elemental data into the array m_phys. More...
 
Array< OneD, int > m_offset_elmt_id
 Array containing the element id m_offset_elmt_id[n] that the n^th consecutive block of data in m_coeffs and m_phys is associated, i.e. for an array of constant expansion size and single shape elements m_phys[n*m_npoints] is the data related to m_exp[m_offset_elmt_id[n]];. More...
 
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
 
BlockMatrixMapShPtr m_blockMat
 
bool m_WaveSpace
 

Detailed Description

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

109  :
110  ExpList(pSession,graph1D)
111  {
112  SetExpType(e1D);
113 
114  int id=0;
117 
118  const SpatialDomains::ExpansionMap &expansions
119  = graph1D->GetExpansions();
120 
121  // For each element in the mesh, create a segment expansion using
122  // the supplied BasisKey and segment geometry.
123  SpatialDomains::ExpansionMap::const_iterator expIt;
124  for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
125  {
126  if ((SegmentGeom = boost::dynamic_pointer_cast<
128  expIt->second->m_geomShPtr)))
129  {
131  ::AllocateSharedPtr(Ba,SegmentGeom);
132  seg->SetElmtId(id++);
133  (*m_exp).push_back(seg);
134  }
135  else
136  {
137  ASSERTL0(false,"dynamic cast to a SegGeom failed");
138  }
139  }
140 
141  // Setup Default optimisation information.
142  int nel = GetExpSize();
145 
146  // Allocate storage for data and populate element offset lists.
148 
151 
154  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
ExpList()
The default constructor.
Definition: ExpList.cpp:93
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:971
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:926
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:909
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1837
boost::shared_ptr< SegExp > SegExpSharedPtr
Definition: SegExp.h:255
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:887
void SetCoeffPhysOffsets(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
Definition: ExpList1D.cpp:685
void CreateCollections(Collections::ImplementationType ImpType=Collections::eNoImpType)
Construct collections of elements containing a single element type and polynomial order from the list...
Definition: ExpList.cpp:2723
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:211
std::map< int, ExpansionShPtr > ExpansionMap
Definition: MeshGraph.h:171
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 177 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().

179  :
180  ExpList(pSession,graph1D)
181  {
182  SetExpType(e1D);
183 
184  int id=0;
187 
188  // Retrieve the list of expansions
189  const SpatialDomains::ExpansionMap &expansions
190  = graph1D->GetExpansions();
191 
192  // Process each expansion in the graph
193  SpatialDomains::ExpansionMap::const_iterator expIt;
194  for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
195  {
196  // Retrieve basis key from expansion
197  LibUtilities::BasisKey bkey = expIt->second->m_basisKeyVector[0];
198 
199  if ((SegmentGeom = boost::dynamic_pointer_cast<
201  expIt->second->m_geomShPtr)))
202  {
204  ::AllocateSharedPtr(bkey, SegmentGeom);
205 
206  // Assign next ID
207  seg->SetElmtId(id++);
208 
209  // Add the expansion
210  (*m_exp).push_back(seg);
211  }
212  else
213  {
214  ASSERTL0(false,"dynamic cast to a SegGeom failed");
215  }
216  }
217 
218  // Setup Default optimisation information.
219  int nel = GetExpSize();
222 
223  // set up offset arrays.
225 
226  if(DeclareCoeffPhysArrays)
227  {
228  // Set up m_coeffs, m_phys.
231  }
232 
235  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
ExpList()
The default constructor.
Definition: ExpList.cpp:93
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:971
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:926
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:909
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1837
boost::shared_ptr< SegExp > SegExpSharedPtr
Definition: SegExp.h:255
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:887
void SetCoeffPhysOffsets(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
Definition: ExpList1D.cpp:685
void CreateCollections(Collections::ImplementationType ImpType=Collections::eNoImpType)
Construct collections of elements containing a single element type and polynomial order from the list...
Definition: ExpList.cpp:2723
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:211
Describes the specification for a Basis.
Definition: Basis.h:50
std::map< int, ExpansionShPtr > ExpansionMap
Definition: MeshGraph.h:171
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 260 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().

265  :
266  ExpList(pSession,graph1D)
267  {
268  int j,id=0;
272  SpatialDomains::CompositeMap::const_iterator compIt;
273 
274  // Retrieve the list of expansions
275  const SpatialDomains::ExpansionMap &expansions
276  = graph1D->GetExpansions(var);
277 
278  // Process each composite region in domain
279  for(compIt = domain.begin(); compIt != domain.end(); ++compIt)
280  {
281  comp = compIt->second;
282 
283  // Process each expansion in the graph
284  for(j = 0; j < compIt->second->size(); ++j)
285  {
286  SpatialDomains::ExpansionMap::const_iterator expIt;
287 
288  if((SegmentGeom = boost::dynamic_pointer_cast<
290  (*compIt->second)[j])))
291  {
292  // Retrieve basis key from expansion and define expansion
293  if((expIt = expansions.find(SegmentGeom->GetGlobalID())) != expansions.end())
294  {
295  LibUtilities::BasisKey bkey = expIt->second->m_basisKeyVector[0];
296 
297  if(SetToOneSpaceDimension)
298  {
299  SpatialDomains::SegGeomSharedPtr OneDSegmentGeom =
300  SegmentGeom->GenerateOneSpaceDimGeom();
301 
303  ::AllocateSharedPtr(bkey, OneDSegmentGeom);
304  }
305  else
306  {
308  ::AllocateSharedPtr(bkey, SegmentGeom);
309  }
310  }
311  else
312  {
313  ASSERTL0(false,"Failed to find basis key");
314  }
315 
316  }
317  else
318  {
319  ASSERTL0(false,"Failed to dynamic cast geometry to SegGeom");
320  }
321 
322 
323  // Assign next ID
324  seg->SetElmtId(id++);
325 
326  // Add the expansion
327  (*m_exp).push_back(seg);
328  }
329  }
330 
331  // Setup Default optimisation information.
332  int nel = GetExpSize();
335 
336  // set up offset arrays.
338 
339  if(DeclareCoeffPhysArrays)
340  {
341  // Set up m_coeffs, m_phys.
344  }
345 
348  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
ExpList()
The default constructor.
Definition: ExpList.cpp:93
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:971
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:926
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:909
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1837
boost::shared_ptr< SegExp > SegExpSharedPtr
Definition: SegExp.h:255
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:887
boost::shared_ptr< GeometryVector > Composite
Definition: MeshGraph.h:111
void SetCoeffPhysOffsets(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
Definition: ExpList1D.cpp:685
void CreateCollections(Collections::ImplementationType ImpType=Collections::eNoImpType)
Construct collections of elements containing a single element type and polynomial order from the list...
Definition: ExpList.cpp:2723
Describes the specification for a Basis.
Definition: Basis.h:50
std::map< int, ExpansionShPtr > ExpansionMap
Definition: MeshGraph.h:171
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 364 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().

368  :
369  ExpList(pSession,graph2D)
370  {
371  SetExpType(e1D);
372 
373  m_graph = graph2D;
374 
375  int j, id=0;
377  SpatialDomains::CompositeMap::const_iterator compIt;
380 
381  // Process each composite region.
382  for(compIt = domain.begin(); compIt != domain.end(); ++compIt)
383  {
384  comp = compIt->second;
385  // Process each expansion in the region.
386  for(j = 0; j < compIt->second->size(); ++j)
387  {
388  if((SegmentGeom = boost::dynamic_pointer_cast<
390  (*compIt->second)[j])))
391  {
392  // Retrieve the basis key from the expansion.
394  = boost::dynamic_pointer_cast<SpatialDomains::MeshGraph2D>(graph2D)->GetEdgeBasisKey(SegmentGeom, variable);
395 
397  ::AllocateSharedPtr(bkey, SegmentGeom);
398 
399  // Add the segment to the expansion list.
400  seg->SetElmtId(id++);
401  (*m_exp).push_back(seg);
402  }
403  else
404  {
405  ASSERTL0(false,"dynamic cast to a SegGeom failed");
406  }
407  }
408 
409  }
410 
411  // Setup Default optimisation information.
412  int nel = GetExpSize();
415 
416  // Allocate storage for data and populate element offset lists.
418 
419  // Set up m_coeffs, m_phys.
420  if(DeclareCoeffPhysArrays)
421  {
424  }
425 
427  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
ExpList()
The default constructor.
Definition: ExpList.cpp:93
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:971
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:926
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:909
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1837
boost::shared_ptr< SegExp > SegExpSharedPtr
Definition: SegExp.h:255
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:887
SpatialDomains::MeshGraphSharedPtr m_graph
Mesh associated with this expansion list.
Definition: ExpList.h:883
boost::shared_ptr< GeometryVector > Composite
Definition: MeshGraph.h:111
void SetCoeffPhysOffsets(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
Definition: ExpList1D.cpp:685
void CreateCollections(Collections::ImplementationType ImpType=Collections::eNoImpType)
Construct collections of elements containing a single element type and polynomial order from the list...
Definition: ExpList.cpp:2723
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:211
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 445 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().

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

710  {
711  }

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

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

Referenced by PostProcess().

826  {
827  int i,j,r;
828 
829  // Get the number of elements in the domain
830  int num_elm = GetExpSize();
831 
832  // initializing the outarray
833  for(i = 0; i < outarray.num_elements(); i++)
834  {
835  outarray[i] = 0.0;
836  }
837 
838  // Make a copy for further modification
839  int x_size = inarray2.num_elements();
840  Array<OneD,NekDouble> x_values_cp(x_size);
841 
842  // Determining the element to which the x belongs
843  Array<OneD,int> x_elm(x_size);
844  for(i = 0; i < x_size; i++ )
845  {
846  x_elm[i] = (int)floor(inarray2[i]/h);
847  }
848 
849  // Clamp indices periodically
850  for(i = 0; i < x_size; i++)
851  {
852  while(x_elm[i] < 0)
853  {
854  x_elm[i] += num_elm;
855  }
856  while(x_elm[i] >= num_elm)
857  {
858  x_elm[i] -= num_elm ;
859  }
860  }
861 
862  // Map the values of x to [-1 1] on its interval
863  for(i = 0; i < x_size; i++)
864  {
865  x_values_cp[i] = (inarray2[i]/h - floor(inarray2[i]/h))*2 - 1.0;
866  }
867 
868  // Evaluate the jocobi polynomials
869  // (Evaluating the base at some points other than the quadrature
870  // points). Should it be added to the base class????
871  Array<TwoD,NekDouble> jacobi_poly(nmodes,x_size);
872  for(i = 0; i < nmodes; i++)
873  {
874  Polylib::jacobfd(x_size,x_values_cp.get(),
875  jacobi_poly.get()+i*x_size,NULL,i,0.0,0.0);
876  }
877 
878  // Evaluate the function values
879  for(r = 0; r < nmodes; r++)
880  {
881  for(j = 0; j < x_size; j++)
882  {
883  int index = ((x_elm[j])*nmodes)+r;
884  outarray[j] += inarray1[index]*jacobi_poly[r][j];
885  }
886  }
887 
888  }
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1837
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 723 of file ExpList1D.cpp.

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

729  {
730  int i,j,r;
731 
732  // get the local element expansion of the elmId element
734 
735  // Get the quadrature points and weights required for integration
736  int quad_npoints = elmExp->GetTotPoints();
737  LibUtilities::PointsKey quadPointsKey(quad_npoints,
738  elmExp->GetPointsType(0));
739  Array<OneD,NekDouble> quad_points
740  = LibUtilities::PointsManager()[quadPointsKey]->GetZ();
741  Array<OneD,NekDouble> quad_weights
742  = LibUtilities::PointsManager()[quadPointsKey]->GetW();
743 
744  // Declare variable for the local kernel breaks
745  int kernel_width = kernel->GetKernelWidth();
746  Array<OneD,NekDouble> local_kernel_breaks(kernel_width+1);
747 
748  // Declare variable for the transformed quadrature points
749  Array<OneD,NekDouble> mapped_quad_points(quad_npoints);
750 
751  // For each evaluation point
752  for(i = 0; i < inarray.num_elements(); i++)
753  {
754  // Move the center of the kernel to the current point
755  kernel->MoveKernelCenter(inarray[i],local_kernel_breaks);
756 
757  // Find the mesh breaks under the kernel support
758  Array<OneD,NekDouble> mesh_breaks;
759  kernel->FindMeshUnderKernel(local_kernel_breaks,h,mesh_breaks);
760 
761  // Sort the total breaks for integration purposes
762  int total_nbreaks = local_kernel_breaks.num_elements() +
763  mesh_breaks.num_elements();
764  // number of the total breaks
765  Array<OneD,NekDouble> total_breaks(total_nbreaks);
766  kernel->Sort(local_kernel_breaks,mesh_breaks,total_breaks);
767 
768  // Integrate the product of kernel and function over the total
769  // breaks
770  NekDouble integral_value = 0.0;
771  for(j = 0; j < total_breaks.num_elements()-1; j++)
772  {
773  NekDouble a = total_breaks[j];
774  NekDouble b = total_breaks[j+1];
775 
776  // Map the quadrature points to the appropriate interval
777  for(r = 0; r < quad_points.num_elements(); r++)
778  {
779  mapped_quad_points[r]
780  = (quad_points[r] + 1.0) * 0.5 * (b - a) + a;
781  }
782 
783  // Evaluate the function at the transformed quadrature
784  // points
785  Array<OneD,NekDouble> u_value(quad_npoints);
786  Array<OneD,NekDouble> coeffs = GetCoeffs();
787 
788  PeriodicEval(coeffs,mapped_quad_points,h,
789  elmExp->GetBasisNumModes(0),u_value);
790 
791  // Evaluate the kernel at the transformed quadrature points
792  Array<OneD,NekDouble> k_value(quad_npoints);
793  kernel->EvaluateKernel(mapped_quad_points,h,k_value);
794 
795  // Integrate
796  for(r = 0; r < quad_npoints; r++)
797  {
798  integral_value += (b - a) * 0.5 * k_value[r]
799  * u_value[r] * quad_weights[r];
800  }
801  }
802  outarray[i] = integral_value/h;
803  }
804  }
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:1775
const boost::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:1858
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:822
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 685 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().

686  {
687  int i;
688 
689  // Set up offset information and array sizes
693 
694  m_ncoeffs = m_npoints = 0;
695 
696  for(i = 0; i < m_exp->size(); ++i)
697  {
699  m_phys_offset [i] = m_npoints;
700  m_offset_elmt_id[i] = i;
701  m_ncoeffs += (*m_exp)[i]->GetNcoeffs();
702  m_npoints += (*m_exp)[i]->GetTotPoints();
703  }
704  }
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:958
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:947
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:961
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:887
Array< OneD, int > m_offset_elmt_id
Array containing the element id m_offset_elmt_id[n] that the n^th consecutive block of data in m_coef...
Definition: ExpList.h:969
void Nektar::MultiRegions::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 1050 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.

1052  {
1053  int i,j,k,e_npoints,offset;
1055  Array<OneD,Array<OneD,NekDouble> > locnormals;
1056  Array<OneD,Array<OneD,NekDouble> > locnormals2;
1058  // Assume whole array is of same coordinate dimension
1059  int coordim = GetCoordim(0);
1060 
1061  ASSERTL1(normals.num_elements() >= coordim,
1062  "Output vector does not have sufficient dimensions to "
1063  "match coordim");
1064 
1065  for (i = 0; i < m_exp->size(); ++i)
1066  {
1068 
1070  loc_exp->GetLeftAdjacentElementExp();
1071 
1072  int edgeNumber = loc_exp->GetLeftAdjacentElementEdge();
1073 
1074  // Get the number of points and normals for this expansion.
1075  e_npoints = (*m_exp)[i]->GetNumPoints(0);
1076 
1077  locnormals = loc_elmt->GetEdgeNormal(edgeNumber);
1078  int e_nmodes = loc_exp->GetBasis(0)->GetNumModes();
1079  int loc_nmodes = loc_elmt->GetBasis(0)->GetNumModes();
1080 
1081  if (e_nmodes != loc_nmodes)
1082  {
1083  if (loc_exp->GetRightAdjacentElementEdge() >= 0)
1084  {
1086  loc_exp->GetRightAdjacentElementExp();
1087 
1088  int EdgeNumber = loc_exp->GetRightAdjacentElementEdge();
1089  // Serial case: right element is connected so we can
1090  // just grab that normal.
1091  locnormals = loc_elmt->GetEdgeNormal(EdgeNumber);
1092 
1093  offset = m_phys_offset[i];
1094 
1095  // Process each point in the expansion.
1096  for (j = 0; j < e_npoints; ++j)
1097  {
1098  // Process each spatial dimension and copy the values
1099  // into the output array.
1100  for (k = 0; k < coordim; ++k)
1101  {
1102  normals[k][offset + j] = -locnormals[k][j];
1103  }
1104  }
1105  }
1106  else
1107  {
1108  // Parallel case: need to interpolate normal.
1109  Array<OneD, Array<OneD, NekDouble> > normal(coordim);
1110 
1111  for (int p = 0; p < coordim; ++p)
1112  {
1113  normal[p] = Array<OneD, NekDouble>(e_npoints,0.0);
1114  LibUtilities::PointsKey to_key =
1115  loc_exp->GetBasis(0)->GetPointsKey();
1116  LibUtilities::PointsKey from_key =
1117  loc_elmt->GetBasis(0)->GetPointsKey();
1118  LibUtilities::Interp1D(from_key,
1119  locnormals[p],
1120  to_key,
1121  normal[p]);
1122  }
1123 
1124  offset = m_phys_offset[i];
1125 
1126  // Process each point in the expansion.
1127  for (j = 0; j < e_npoints; ++j)
1128  {
1129  // Process each spatial dimension and copy the values
1130  // into the output array.
1131  for (k = 0; k < coordim; ++k)
1132  {
1133  normals[k][offset + j] = normal[k][j];
1134  }
1135  }
1136  }
1137  }
1138  else
1139  {
1140  // Get the physical data offset for this expansion.
1141  offset = m_phys_offset[i];
1142 
1143  // Process each point in the expansion.
1144  for (j = 0; j < e_npoints; ++j)
1145  {
1146  // Process each spatial dimension and copy the values
1147  // into the output array.
1148  for (k = 0; k < coordim; ++k)
1149  {
1150  normals[k][offset + j] = locnormals[k][j];
1151  }
1152  }
1153  }
1154  }
1155  }
Expansion2DSharedPtr GetLeftAdjacentElementExp() const
Definition: Expansion1D.h:116
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:947
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:961
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:1735
#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 1160 of file ExpList1D.cpp.

1161  {
1162 // Array<OneD, int> NumShape(1,0);
1163 // NumShape[0] = GetExpSize();
1164 //
1165 // int one = 1;
1166 // m_globalOptParam = MemoryManager<NekOptimize::GlobalOptParam>
1167 // ::AllocateSharedPtr(m_session,one,NumShape);
1168  }
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 925 of file ExpList1D.cpp.

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

926  {
927  int i, j;
928  for (i = 0; i < m_exp->size(); ++i)
929  {
930  for (j = 0; j < (*m_exp)[i]->GetNverts(); ++j)
931  {
932  (*m_exp)[i]->ComputeVertexNormal(j);
933  }
934  }
935  }
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:947
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 947 of file ExpList1D.cpp.

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

952  {
953  int i,j,k,e_npoints,offset;
954  Array<OneD,NekDouble> normals;
955  NekDouble Vn;
956 
957  // Assume whole array is of same coordimate dimension
958  int coordim = GetCoordim(0);
959 
960  ASSERTL1(Vec.num_elements() >= coordim,
961  "Input vector does not have sufficient dimensions to "
962  "match coordim");
963 
964  // Process each expansion
965  for(i = 0; i < m_exp->size(); ++i)
966  {
967  // Get the number of points in the expansion and the normals.
968  e_npoints = (*m_exp)[i]->GetNumPoints(0);
969  normals = (*m_exp)[i]->GetPhysNormals();
970 
971  // Get the physical data offset of the expansion in m_phys.
972  offset = m_phys_offset[i];
973 
974  // Compute each data point.
975  for(j = 0; j < e_npoints; ++j)
976  {
977  // Calculate normal velocity.
978  Vn = 0.0;
979  for(k = 0; k < coordim; ++k)
980  {
981  Vn += Vec[k][offset+j]*normals[k*e_npoints + j];
982  }
983 
984  // Upwind based on direction of normal velocity.
985  if(Vn > 0.0)
986  {
987  Upwind[offset + j] = Fwd[offset + j];
988  }
989  else
990  {
991  Upwind[offset + j] = Bwd[offset + j];
992  }
993  }
994  }
995  }
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:947
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:961
double NekDouble
int GetCoordim(int eid)
This function returns the dimension of the coordinates of the element eid.
Definition: ExpList.h:1735
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
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 1010 of file ExpList1D.cpp.

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

1015  {
1016  int i,j,e_npoints,offset;
1017  Array<OneD,NekDouble> normals;
1018 
1019  // Process each expansion.
1020  for(i = 0; i < m_exp->size(); ++i)
1021  {
1022  // Get the number of points and the data offset.
1023  e_npoints = (*m_exp)[i]->GetNumPoints(0);
1024  offset = m_phys_offset[i];
1025 
1026  // Process each point in the expansion.
1027  for(j = 0; j < e_npoints; ++j)
1028  {
1029  // Upwind based on one-dimensional velocity.
1030  if(Vn[offset + j] > 0.0)
1031  {
1032  Upwind[offset + j] = Fwd[offset + j];
1033  }
1034  else
1035  {
1036  Upwind[offset + j] = Bwd[offset + j];
1037  }
1038  }
1039  }
1040  }
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:947
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:961
void Nektar::MultiRegions::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 1179 of file ExpList1D.cpp.

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

Member Data Documentation

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

Definition at line 185 of file ExpList1D.h.

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

Definition at line 187 of file ExpList1D.h.