Nektar++
Public Member Functions | Protected Member Functions | Private Member Functions | List of all members
Nektar::MultiRegions::ExpList3DHomogeneous1D Class Reference

Abstraction of a two-dimensional multi-elemental expansion which is merely a collection of local expansions. More...

#include <ExpList3DHomogeneous1D.h>

Inheritance diagram for Nektar::MultiRegions::ExpList3DHomogeneous1D:
[legend]

Public Member Functions

 ExpList3DHomogeneous1D ()
 Default constructor. More...
 
 ExpList3DHomogeneous1D (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::BasisKey &HomoBasis, const NekDouble lhom, const bool useFFT, const bool dealiasing)
 
 ExpList3DHomogeneous1D (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::BasisKey &HomoBasis, const NekDouble lhom, const bool useFFT, const bool dealiasing, const SpatialDomains::MeshGraphSharedPtr &graph2D, const std::string &var="DefaultVar", const Collections::ImplementationType ImpType=Collections::eNoImpType)
 Sets up a list of local expansions based on an input mesh. More...
 
 ExpList3DHomogeneous1D (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::BasisKey &HomoBasis, const NekDouble lhom, const bool useFFT, const bool dealiasing, const SpatialDomains::ExpansionInfoMap &expansions, const Collections::ImplementationType ImpType=Collections::eNoImpType)
 Sets up a list of local expansions based on an mesh expansion. More...
 
 ExpList3DHomogeneous1D (const ExpList3DHomogeneous1D &In, const bool DeclarePlanesSetCoeffPhys=true)
 Copy constructor. More...
 
 ExpList3DHomogeneous1D (const ExpList3DHomogeneous1D &In, const std::vector< unsigned int > &eIDs, const bool DeclarePlanesSetCoeffPhys=true, const Collections::ImplementationType ImpType=Collections::eNoImpType)
 Constructor copying only elements defined in eIds. More...
 
 ~ExpList3DHomogeneous1D () override
 Destructor. More...
 
- Public Member Functions inherited from Nektar::MultiRegions::ExpListHomogeneous1D
 ExpListHomogeneous1D (const ExpansionType type)
 Default constructor. More...
 
 ExpListHomogeneous1D (const ExpansionType type, const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::BasisKey &HomoBasis, const NekDouble lz, const bool useFFT, const bool dealiasing)
 
 ExpListHomogeneous1D (const ExpListHomogeneous1D &In)
 Copy constructor. More...
 
 ExpListHomogeneous1D (const ExpListHomogeneous1D &In, const std::vector< unsigned int > &eIDs, const Collections::ImplementationType ImpType=Collections::eNoImpType)
 
 ~ExpListHomogeneous1D () override
 Destructor. More...
 
void Homogeneous1DTrans (const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool IsForwards, bool Shuff=true, bool UnShuff=true)
 
LibUtilities::BasisSharedPtr GetHomogeneousBasis (void)
 
void PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
 
void PhysDeriv (Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
 
ExpListSharedPtrGetPlane (int n)
 
- Public Member Functions inherited from Nektar::MultiRegions::ExpList
 ExpList (const ExpansionType Type=eNoType)
 The default constructor using a type. More...
 
 ExpList (const ExpList &in, const bool DeclareCoeffPhysArrays=true)
 The copy constructor. More...
 
 ExpList (const ExpList &in, const std::vector< unsigned int > &eIDs, const bool DeclareCoeffPhysArrays=true, const Collections::ImplementationType ImpType=Collections::eNoImpType)
 Constructor copying only elements defined in eIds. More...
 
 ExpList (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &graph, const bool DeclareCoeffPhysArrays=true, const std::string &var="DefaultVar", const Collections::ImplementationType ImpType=Collections::eNoImpType)
 Generate an ExpList from a meshgraph graph and session file. More...
 
 ExpList (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::ExpansionInfoMap &expansions, const bool DeclareCoeffPhysArrays=true, const Collections::ImplementationType ImpType=Collections::eNoImpType)
 Sets up a list of local expansions based on an expansion Map. More...
 
 ExpList (const SpatialDomains::PointGeomSharedPtr &geom)
 Specialised constructors for 0D Expansions Wrapper around LocalRegion::PointExp - used in PrePacing.cpp. More...
 
 ExpList (const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, const ExpListSharedPtr > &bndConstraint, const Array< OneD, const SpatialDomains ::BoundaryConditionShPtr > &bndCond, const LocalRegions::ExpansionVector &locexp, const SpatialDomains::MeshGraphSharedPtr &graph, const LibUtilities::CommSharedPtr &comm, const bool DeclareCoeffPhysArrays=true, const std::string variable="DefaultVar", const Collections::ImplementationType ImpType=Collections::eNoImpType)
 Generate expansions for the trace space expansions used in DisContField. More...
 
 ExpList (const LibUtilities::SessionReaderSharedPtr &pSession, const LocalRegions::ExpansionVector &locexp, const SpatialDomains::MeshGraphSharedPtr &graph, const bool DeclareCoeffPhysArrays, const std::string variable, const Collections::ImplementationType ImpType=Collections::eNoImpType)
 Generate an trace ExpList from a meshgraph graph and session file. More...
 
 ExpList (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::CompositeMap &domain, const SpatialDomains::MeshGraphSharedPtr &graph, const bool DeclareCoeffPhysArrays=true, const std::string variable="DefaultVar", bool SetToOneSpaceDimension=false, const LibUtilities::CommSharedPtr comm=LibUtilities::CommSharedPtr(), const Collections::ImplementationType ImpType=Collections::eNoImpType)
 Constructor based on domain information only for 1D & 2D boundary conditions. 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...
 
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 (int i, NekDouble val)
 Set the i th value of m_phys to value val. 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...
 
void MultiplyByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 multiply the metric jacobi and quadrature weights More...
 
void DivideByQuadratureMetric (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Divided by the metric jacobi and quadrature weights. More...
 
void IProductWRTBase (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 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 IProductWRTDirectionalDerivBase (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 Directional derivative along a given 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 of all local expansion modes \(\phi_n^e(\boldsymbol{x})\). More...
 
void FwdTransLocalElmt (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)
 
void ExponentialFilter (Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
 
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)
 
void MultiplyByMassMatrix (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void SmoothField (Array< OneD, NekDouble > &field)
 Smooth a field across elements. More...
 
GlobalLinSysKey HelmSolve (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff=StdRegions::NullVarCoeffMap, const MultiRegions::VarFactorsMap &varfactors=MultiRegions::NullVarFactorsMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray, const bool PhysSpaceForcing=true)
 Solve helmholtz problem. More...
 
GlobalLinSysKey LinearAdvectionDiffusionReactionSolve (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff=StdRegions::NullVarCoeffMap, const MultiRegions::VarFactorsMap &varfactors=MultiRegions::NullVarFactorsMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray, const bool PhysSpaceForcing=true)
 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, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 Solve Advection Diffusion Reaction. More...
 
void FwdTransBndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void BwdTrans (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 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 GetCoords (const int eid, Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1=NullNekDouble1DArray, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
 
void HomogeneousFwdTrans (const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
 
void HomogeneousBwdTrans (const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
 
void DealiasedProd (const int num_dofs, const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray)
 
void DealiasedDotProd (const int num_dofs, const Array< OneD, Array< OneD, NekDouble > > &inarray1, const Array< OneD, Array< OneD, NekDouble > > &inarray2, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
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...
 
int GetShapeDimension ()
 This function returns the dimension of the shape of the element eid. 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 (const Array< OneD, NekDouble > coeffs)
 Fill Bnd Condition expansion from the values stored in expansion. More...
 
void FillBndCondFromField (const int nreg, const Array< OneD, NekDouble > coeffs)
 Fill Bnd Condition expansion in nreg from the values stored in expansion. More...
 
void LocalToGlobal (bool useComm=true)
 Gathers the global coefficients \(\boldsymbol{\hat{u}}_g\) from the local coefficients \(\boldsymbol{\hat{u}}_l\). More...
 
void LocalToGlobal (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool useComm=true)
 
void GlobalToLocal (void)
 Scatters from the global coefficients \(\boldsymbol{\hat{u}}_g\) to the local coefficients \(\boldsymbol{\hat{u}}_l\). More...
 
void GlobalToLocal (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
NekDouble GetCoeff (int i)
 Get the i th value (coefficient) of m_coeffs. More...
 
NekDouble GetCoeffs (int i)
 Get the i th value (coefficient) of m_coeffs. More...
 
const Array< OneD, const NekDouble > & GetPhys () const
 This function returns (a reference to) the array \(\boldsymbol{u}_l\) (implemented as m_phys) containing the function \(u^{\delta}(\boldsymbol{x})\) evaluated at the quadrature points. More...
 
NekDouble Linf (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
 This function calculates the \(L_\infty\) error of the global spectral/hp element approximation. More...
 
NekDouble L2 (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
 This function calculates the \(L_\infty\) error of the global 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 ()
 Calculates the \(H^1\) error of the global spectral/hp element approximation. More...
 
NekDouble Integral (const Array< OneD, const NekDouble > &inarray)
 
NekDouble VectorFlux (const Array< OneD, Array< OneD, 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 associated with the homogeneous expansion. More...
 
NekDouble GetHomoLen (void)
 This function returns the Width of homogeneous direction associated with the homogeneous expansion. More...
 
void SetHomoLen (const NekDouble lhom)
 This function sets the Width of homogeneous direction associated 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...
 
size_t GetNumElmts (void)
 This function returns the number of elements in the expansion which may be different for a homogeoenous extended expansionp. More...
 
const std::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::ExpansionSharedPtrGetExpFromGeomId (int n)
 This function returns (a shared pointer to) the local elemental expansion of the \(n^{\mathrm{th}}\) element given a global geometry ID. 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 cachedId=-1, NekDouble maxDistance=1e6)
 This function returns the index of the local elemental expansion containing the arbitrary point given by gloCoord, within a distance tolerance of tol. More...
 
int GetExpIndex (const Array< OneD, const NekDouble > &gloCoords, Array< OneD, NekDouble > &locCoords, NekDouble tol=0.0, bool returnNearestElmt=false, int cachedId=-1, NekDouble maxDistance=1e6)
 
NekDouble PhysEvaluate (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &phys)
 
int GetCoeff_Offset (int n) const
 Get the start offset position for a local contiguous list of coeffs correspoinding to element n. More...
 
int GetPhys_Offset (int n) const
 Get the start offset position for a local contiguous list of quadrature points in a full array correspoinding to element n. More...
 
Array< OneD, NekDouble > & UpdateCoeffs ()
 This function returns (a reference to) the array \(\boldsymbol{\hat{u}}_l\) (implemented as m_coeffs) containing all local expansion coefficients. More...
 
Array< OneD, NekDouble > & UpdatePhys ()
 This function returns (a reference to) the array \(\boldsymbol{u}_l\) (implemented as m_phys) containing the function \(u^{\delta}(\boldsymbol{x})\) evaluated at the quadrature points. More...
 
void PhysDeriv (Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
 
void PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 This function discretely evaluates the derivative of a function \(f(\boldsymbol{x})\) on the domain consisting of all elements of the expansion. More...
 
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
 
void Curl (Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
 
void CurlCurl (Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
 
void PhysDirectionalDeriv (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
void GetMovingFrames (const SpatialDomains::GeomMMF MMFdir, const Array< OneD, const NekDouble > &CircCentre, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
const Array< OneD, const std::shared_ptr< ExpList > > & GetBndCondExpansions ()
 
const Array< OneD, const NekDouble > & GetBndCondBwdWeight ()
 Get the weight value for boundary conditions. More...
 
void SetBndCondBwdWeight (const int index, const NekDouble value)
 Set the weight value for boundary conditions. More...
 
std::shared_ptr< ExpList > & UpdateBndCondExpansion (int i)
 
void Upwind (const Array< OneD, const NekDouble > &Vn, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
 
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)
 
std::shared_ptr< ExpList > & GetTrace ()
 
std::shared_ptr< AssemblyMapDG > & GetTraceMap (void)
 
const Array< OneD, const int > & GetTraceBndMap (void)
 
void GetNormals (Array< OneD, Array< OneD, NekDouble > > &normals)
 
void GetElmtNormalLength (Array< OneD, NekDouble > &lengthsFwd, Array< OneD, NekDouble > &lengthsBwd)
 Get the length of elements in boundary normal direction. More...
 
void GetBwdWeight (Array< OneD, NekDouble > &weightAver, Array< OneD, NekDouble > &weightJump)
 Get the weight value for boundary conditions for boundary average and jump calculations. More...
 
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, bool FillBnd=true, bool PutFwdInBwdOnBCs=false, bool DoExchange=true)
 
void FillBwdWithBoundCond (const Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd, bool PutFwdInBwdOnBCs=false)
 
void AddTraceQuadPhysToField (const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field)
 Add Fwd and Bwd value to field, a reverse procedure of GetFwdBwdTracePhys. More...
 
void AddTraceQuadPhysToOffDiag (const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field)
 
void GetLocTraceFromTracePts (const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &locTraceFwd, Array< OneD, NekDouble > &locTraceBwd)
 
void FillBwdWithBwdWeight (Array< OneD, NekDouble > &weightave, Array< OneD, NekDouble > &weightjmp)
 Fill Bwd with boundary conditions. More...
 
void PeriodicBwdCopy (const Array< OneD, const NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
 Copy and fill the Periodic boundaries. More...
 
const std::vector< bool > & GetLeftAdjacentFaces (void) const
 
void ExtractTracePhys (Array< OneD, NekDouble > &outarray)
 
void ExtractTracePhys (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & GetBndConditions ()
 
Array< OneD, SpatialDomains::BoundaryConditionShPtr > & UpdateBndConditions ()
 
void EvaluateBoundaryConditions (const NekDouble time=0.0, const std::string varName="", const NekDouble=NekConstants::kNekUnsetDouble, const NekDouble=NekConstants::kNekUnsetDouble)
 
void GeneralMatrixOp (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function calculates the result of the multiplication of a matrix of type specified by mkey with a vector given by inarray. More...
 
void SetUpPhysNormals ()
 
void GetBoundaryToElmtMap (Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
 
virtual void GetBndElmtExpansion (int i, std::shared_ptr< ExpList > &result, const bool DeclareCoeffPhysArrays=true)
 
void ExtractElmtToBndPhys (int i, const Array< OneD, NekDouble > &elmt, Array< OneD, NekDouble > &boundary)
 
void ExtractPhysToBndElmt (int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bndElmt)
 
void ExtractPhysToBnd (int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bnd)
 
void GetBoundaryNormals (int i, Array< OneD, Array< OneD, NekDouble > > &normals)
 
void GeneralGetFieldDefinitions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef, int NumHomoDir=0, Array< OneD, LibUtilities::BasisSharedPtr > &HomoBasis=LibUtilities::NullBasisSharedPtr1DArray, std::vector< NekDouble > &HomoLen=LibUtilities::NullNekDoubleVector, bool homoStrips=false, std::vector< unsigned int > &HomoSIDs=LibUtilities::NullUnsignedIntVector, std::vector< unsigned int > &HomoZIDs=LibUtilities::NullUnsignedIntVector, std::vector< unsigned int > &HomoYIDs=LibUtilities::NullUnsignedIntVector)
 
std::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 std::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, std::unordered_map< int, int > zIdToPlane=std::unordered_map< int, int >())
 Extract the data in fielddata into the coeffs. More...
 
void ExtractCoeffsFromFile (const std::string &fileName, LibUtilities::CommSharedPtr comm, const std::string &varName, Array< OneD, NekDouble > &coeffs)
 
void GenerateElementVector (const int ElementID, const NekDouble scalar1, const NekDouble scalar2, Array< OneD, NekDouble > &outarray)
 Generate vector v such that v[i] = scalar1 if i is in the element < ElementID. Otherwise, v[i] = scalar2. More...
 
std::shared_ptr< ExpListGetSharedThisPtr ()
 Returns a shared pointer to the current object. More...
 
std::shared_ptr< LibUtilities::SessionReaderGetSession () const
 Returns the session object. More...
 
std::shared_ptr< LibUtilities::CommGetComm () const
 Returns the comm object. More...
 
SpatialDomains::MeshGraphSharedPtr GetGraph ()
 
LibUtilities::BasisSharedPtr GetHomogeneousBasis (void)
 
std::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)
 
int GetPoolCount (std::string)
 
void UnsetGlobalLinSys (GlobalLinSysKey, bool)
 
LibUtilities::NekManager< GlobalLinSysKey, GlobalLinSys > & GetGlobalLinSysManager (void)
 
const Array< OneD, const std::pair< int, int > > & GetCoeffsToElmt () const
 Get m_coeffs to elemental value map. More...
 
void AddTraceJacToElmtJac (const Array< OneD, const DNekMatSharedPtr > &FwdMat, const Array< OneD, const DNekMatSharedPtr > &BwdMat, Array< OneD, DNekMatSharedPtr > &fieldMat)
 inverse process of v_GetFwdBwdTracePhys. Given Trace integration of Fwd and Bwd Jacobian, with dimension NtotalTrace*TraceCoef*TracePhys. return Elemental Jacobian matrix with dimension NtotalElement*ElementCoef*ElementPhys. More...
 
void GetMatIpwrtDeriveBase (const Array< OneD, const Array< OneD, NekDouble > > &inarray, const int nDirctn, Array< OneD, DNekMatSharedPtr > &mtxPerVar)
 
void GetMatIpwrtDeriveBase (const TensorOfArray3D< NekDouble > &inarray, Array< OneD, DNekMatSharedPtr > &mtxPerVar)
 
void GetDiagMatIpwrtBase (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, DNekMatSharedPtr > &mtxPerVar)
 
void AddTraceIntegralToOffDiag (const Array< OneD, const NekDouble > &FwdFlux, const Array< OneD, const NekDouble > &BwdFlux, Array< OneD, NekDouble > &outarray)
 
void AddRightIPTPhysDerivBase (const int dir, const Array< OneD, const DNekMatSharedPtr > ElmtJacQuad, Array< OneD, DNekMatSharedPtr > ElmtJacCoef)
 
void AddRightIPTBaseMatrix (const Array< OneD, const DNekMatSharedPtr > ElmtJacQuad, Array< OneD, DNekMatSharedPtr > ElmtJacCoef)
 
const LocTraceToTraceMapSharedPtrGetLocTraceToTraceMap () const
 
std::vector< bool > & GetLeftAdjacentTraces (void)
 
const std::unordered_map< int, int > & GetElmtToExpId (void)
 This function returns the map of index inside m_exp to geom id. More...
 
int GetElmtToExpId (int elmtId)
 This function returns the index inside m_exp for a given geom id. More...
 

Protected Member Functions

void SetCoeffPhys (void)
 Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m_coeff and m_phys. More...
 
void v_GetCoords (Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1, Array< OneD, NekDouble > &coord_2) override
 
void v_GetCoords (const int eid, Array< OneD, NekDouble > &xc0, Array< OneD, NekDouble > &xc1, Array< OneD, NekDouble > &xc2) final
 
void v_WriteTecplotConnectivity (std::ostream &outfile, int expansion) override
 
void v_WriteVtkPieceHeader (std::ostream &outfile, int expansion, int istrip) override
 
NekDouble v_L2 (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray) override
 
Array< OneD, const NekDoublev_HomogeneousEnergy (void) override
 
void v_GetPeriodicEntities (PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces) override
 
- Protected Member Functions inherited from Nektar::MultiRegions::ExpListHomogeneous1D
DNekBlkMatSharedPtr GenHomogeneous1DBlockMatrix (Homogeneous1DMatType mattype) const
 
DNekBlkMatSharedPtr GetHomogeneous1DBlockMatrix (Homogeneous1DMatType mattype) const
 
NekDouble GetSpecVanVisc (const int k)
 
void v_SetHomo1DSpecVanVisc (Array< OneD, NekDouble > visc) override
 
size_t v_GetNumElmts (void) override
 
LibUtilities::BasisSharedPtr v_GetHomogeneousBasis (void) override
 
void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_FwdTransLocalElmt (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_FwdTransBndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_IProductWRTDerivBase (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &outarray) override
 
std::vector< LibUtilities::FieldDefinitionsSharedPtrv_GetFieldDefinitions (void) override
 
void v_GetFieldDefinitions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef) override
 
void v_AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata) override
 
void v_AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, Array< OneD, NekDouble > &coeffs) override
 
void v_ExtractDataToCoeffs (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs, std::unordered_map< int, int > zIdToPlane) override
 Extract data from raw field data into expansion list. More...
 
void v_ExtractCoeffsToCoeffs (const std::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs) override
 
void v_WriteVtkPieceData (std::ostream &outfile, int expansion, std::string var) override
 
void v_PhysInterp1DScaled (const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_PhysGalerkinProjection1DScaled (const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
 
void v_HomogeneousFwdTrans (const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true) override
 
void v_HomogeneousBwdTrans (const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true) override
 
void v_DealiasedProd (const int num_dofs, const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray) override
 
void v_DealiasedDotProd (const int num_dofs, const Array< OneD, Array< OneD, NekDouble > > &inarray1, const Array< OneD, Array< OneD, NekDouble > > &inarray2, Array< OneD, Array< OneD, NekDouble > > &outarray) override
 
void v_PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2) override
 
void v_PhysDeriv (Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d) override
 
LibUtilities::TranspositionSharedPtr v_GetTransposition (void) override
 
Array< OneD, const unsigned int > v_GetZIDs (void) override
 
ExpListSharedPtrv_GetPlane (int n) override
 
NekDouble v_GetHomoLen (void) override
 
void v_SetHomoLen (const NekDouble lhom) override
 
NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray) override
 
- Protected Member Functions inherited from Nektar::MultiRegions::ExpList
std::shared_ptr< DNekMatGenGlobalMatrixFull (const GlobalLinSysKey &mkey, const std::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)
 
std::shared_ptr< GlobalMatrixGenGlobalMatrix (const GlobalMatrixKey &mkey, const std::shared_ptr< AssemblyMapCG > &locToGloMap)
 Generates a global matrix from the given key and map. More...
 
void GlobalEigenSystem (const std::shared_ptr< DNekMat > &Gmat, Array< OneD, NekDouble > &EigValsReal, Array< OneD, NekDouble > &EigValsImag, Array< OneD, NekDouble > &EigVecs=NullNekDouble1DArray)
 
std::shared_ptr< GlobalLinSysGenGlobalLinSys (const GlobalLinSysKey &mkey, const std::shared_ptr< AssemblyMapCG > &locToGloMap)
 This operation constructs the global linear system of type mkey. More...
 
std::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...
 
virtual size_t v_GetNumElmts (void)
 
virtual const Array< OneD, const std::shared_ptr< ExpList > > & v_GetBndCondExpansions (void)
 
virtual const Array< OneD, const NekDouble > & v_GetBndCondBwdWeight ()
 
virtual void v_SetBndCondBwdWeight (const int index, const NekDouble value)
 
virtual std::shared_ptr< ExpList > & v_UpdateBndCondExpansion (int i)
 
virtual void v_Upwind (const Array< OneD, const Array< OneD, NekDouble > > &Vec, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
 
virtual 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)
 
virtual std::shared_ptr< ExpList > & v_GetTrace ()
 
virtual std::shared_ptr< AssemblyMapDG > & v_GetTraceMap ()
 
virtual const Array< OneD, const int > & v_GetTraceBndMap ()
 
virtual const std::shared_ptr< LocTraceToTraceMap > & v_GetLocTraceToTraceMap (void) const
 
virtual std::vector< bool > & v_GetLeftAdjacentTraces (void)
 
virtual void v_GetNormals (Array< OneD, Array< OneD, NekDouble > > &normals)
 Populate normals with the normals of all expansions. More...
 
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, bool FillBnd=true, bool PutFwdInBwdOnBCs=false, bool DoExchange=true)
 
virtual void v_FillBwdWithBoundCond (const Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd, bool PutFwdInBwdOnBCs)
 
virtual void v_AddTraceQuadPhysToField (const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field)
 
virtual void v_AddTraceQuadPhysToOffDiag (const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field)
 
virtual void v_GetLocTraceFromTracePts (const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &locTraceFwd, Array< OneD, NekDouble > &locTraceBwd)
 
virtual void v_FillBwdWithBwdWeight (Array< OneD, NekDouble > &weightave, Array< OneD, NekDouble > &weightjmp)
 
virtual void v_PeriodicBwdCopy (const Array< OneD, const NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
 
virtual const std::vector< bool > & v_GetLeftAdjacentFaces (void) const
 
virtual void v_ExtractTracePhys (Array< OneD, NekDouble > &outarray)
 
virtual void v_ExtractTracePhys (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_MultiplyByInvMassMatrix (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual GlobalLinSysKey v_HelmSolve (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff, const MultiRegions::VarFactorsMap &varfactors, const Array< OneD, const NekDouble > &dirForcing, const bool PhysSpaceForcing)
 
virtual GlobalLinSysKey v_LinearAdvectionDiffusionReactionSolve (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff, const MultiRegions::VarFactorsMap &varfactors, const Array< OneD, const NekDouble > &dirForcing, const bool PhysSpaceForcing)
 
virtual void v_LinearAdvectionReactionSolve (const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 
virtual void v_ImposeDirichletConditions (Array< OneD, NekDouble > &outarray)
 
virtual void v_FillBndCondFromField (const Array< OneD, NekDouble > coeffs)
 
virtual void v_FillBndCondFromField (const int nreg, const Array< OneD, NekDouble > coeffs)
 
virtual void v_Reset ()
 Reset geometry information, metrics, matrix managers and geometry information. More...
 
virtual void v_LocalToGlobal (bool UseComm)
 
virtual void v_LocalToGlobal (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool UseComm)
 
virtual void v_GlobalToLocal (void)
 
virtual void v_GlobalToLocal (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_FwdTransLocalElmt (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_FwdTransBndConstrained (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)
 
virtual void v_IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTDerivBase (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_GetCoords (Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1=NullNekDouble1DArray, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
 
virtual void v_GetCoords (const int eid, Array< OneD, NekDouble > &xc0, Array< OneD, NekDouble > &xc1, Array< OneD, NekDouble > &xc2)
 
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_Curl (Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
 
virtual void v_CurlCurl (Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
 
virtual void v_PhysDirectionalDeriv (const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_GetMovingFrames (const SpatialDomains::GeomMMF MMFdir, const Array< OneD, const NekDouble > &CircCentre, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual void v_HomogeneousFwdTrans (const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
 
virtual void v_HomogeneousBwdTrans (const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
 
virtual void v_DealiasedProd (const int num_dofs, const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray)
 
virtual void v_DealiasedDotProd (const int num_dofs, const Array< OneD, Array< OneD, NekDouble > > &inarray1, const Array< OneD, Array< OneD, NekDouble > > &inarray2, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
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_SetUpPhysNormals ()
 : Set up a normal along the trace elements between two elements at elemental level More...
 
virtual void v_GetBoundaryToElmtMap (Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
 
virtual void v_GetBndElmtExpansion (int i, std::shared_ptr< ExpList > &result, const bool DeclareCoeffPhysArrays)
 
virtual void v_ExtractElmtToBndPhys (const int i, const Array< OneD, NekDouble > &elmt, Array< OneD, NekDouble > &boundary)
 
virtual void v_ExtractPhysToBndElmt (const int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bndElmt)
 
virtual void v_ExtractPhysToBnd (const int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bnd)
 
virtual void v_GetBoundaryNormals (int i, Array< OneD, Array< OneD, NekDouble > > &normals)
 
virtual std::vector< LibUtilities::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, std::unordered_map< int, int > zIdToPlane)
 Extract data from raw field data into expansion list. More...
 
virtual void v_ExtractCoeffsToCoeffs (const std::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 void v_WriteVtkPieceHeader (std::ostream &outfile, int expansion, int istrip)
 
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 NekDouble v_VectorFlux (const Array< OneD, Array< OneD, NekDouble > > &inarray)
 
virtual Array< OneD, const NekDoublev_HomogeneousEnergy (void)
 
virtual LibUtilities::TranspositionSharedPtr v_GetTransposition (void)
 
virtual NekDouble v_GetHomoLen (void)
 
virtual void v_SetHomoLen (const NekDouble lhom)
 
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)
 
virtual int v_GetPoolCount (std::string)
 
virtual void v_UnsetGlobalLinSys (GlobalLinSysKey, bool)
 
virtual LibUtilities::NekManager< GlobalLinSysKey, GlobalLinSys > & v_GetGlobalLinSysManager (void)
 
void ExtractFileBCs (const std::string &fileName, LibUtilities::CommSharedPtr comm, const std::string &varName, const std::shared_ptr< ExpList > locExpList)
 
virtual const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & v_GetBndConditions ()
 
virtual Array< OneD, SpatialDomains::BoundaryConditionShPtr > & v_UpdateBndConditions ()
 
virtual void v_EvaluateBoundaryConditions (const NekDouble time=0.0, const std::string varName="", const NekDouble x2_in=NekConstants::kNekUnsetDouble, const NekDouble x3_in=NekConstants::kNekUnsetDouble)
 
virtual std::map< int, RobinBCInfoSharedPtrv_GetRobinBCInfo (void)
 
virtual void v_GetPeriodicEntities (PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces)
 
virtual LibUtilities::BasisSharedPtr v_GetHomogeneousBasis (void)
 
virtual void v_SetHomo1DSpecVanVisc (Array< OneD, NekDouble > visc)
 
virtual std::shared_ptr< ExpList > & v_GetPlane (int n)
 
virtual void v_AddTraceIntegralToOffDiag (const Array< OneD, const NekDouble > &FwdFlux, const Array< OneD, const NekDouble > &BwdFlux, Array< OneD, NekDouble > &outarray)
 

Private Member Functions

void GenExpList3DHomogeneous1D (const SpatialDomains::ExpansionInfoMap &expansions, const Collections::ImplementationType ImpType)
 

Additional Inherited Members

- Public Attributes inherited from Nektar::MultiRegions::ExpListHomogeneous1D
LibUtilities::TranspositionSharedPtr m_transposition
 
LibUtilities::CommSharedPtr m_StripZcomm
 
- 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::ExpListHomogeneous1D
bool m_useFFT
 FFT variables. More...
 
LibUtilities::NektarFFTSharedPtr m_FFT
 
LibUtilities::NektarFFTSharedPtr m_FFT_deal
 
Array< OneD, NekDoublem_tmpIN
 
Array< OneD, NekDoublem_tmpOUT
 
LibUtilities::BasisSharedPtr m_homogeneousBasis
 Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m_coeff and m_phys. More...
 
NekDouble m_lhom
 Width of homogeneous direction. More...
 
Homo1DBlockMatrixMapShPtr m_homogeneous1DBlockMat
 
Array< OneD, ExpListSharedPtrm_planes
 
- Protected Attributes inherited from Nektar::MultiRegions::ExpList
ExpansionType m_expType
 Exapnsion type. More...
 
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...
 
std::shared_ptr< LocalRegions::ExpansionVectorm_exp
 The list of local expansions. More...
 
Collections::CollectionVector m_collections
 
std::vector< bool > m_collectionsDoInit
 Vector of bools to act as an initialise on first call flag. More...
 
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, std::pair< int, int > > m_coeffsToElmt
 m_coeffs to elemental value map More...
 
BlockMatrixMapShPtr m_blockMat
 
bool m_WaveSpace
 
std::unordered_map< int, int > m_elmtToExpId
 Mapping from geometry ID of element to index inside m_exp. More...
 

Detailed Description

Abstraction of a two-dimensional multi-elemental expansion which is merely a collection of local expansions.

Definition at line 57 of file ExpList3DHomogeneous1D.h.

Constructor & Destructor Documentation

◆ ExpList3DHomogeneous1D() [1/6]

Nektar::MultiRegions::ExpList3DHomogeneous1D::ExpList3DHomogeneous1D ( )

Default constructor.

Definition at line 44 of file ExpList3DHomogeneous1D.cpp.

45{
46}
ExpListHomogeneous1D(const ExpansionType type)
Default constructor.

◆ ExpList3DHomogeneous1D() [2/6]

Nektar::MultiRegions::ExpList3DHomogeneous1D::ExpList3DHomogeneous1D ( const LibUtilities::SessionReaderSharedPtr pSession,
const LibUtilities::BasisKey HomoBasis,
const NekDouble  lhom,
const bool  useFFT,
const bool  dealiasing 
)

Definition at line 48 of file ExpList3DHomogeneous1D.cpp.

52 : ExpListHomogeneous1D(e3DH1D, pSession, HomoBasis, lhom, useFFT,
53 dealiasing)
54{
55}

◆ ExpList3DHomogeneous1D() [3/6]

Nektar::MultiRegions::ExpList3DHomogeneous1D::ExpList3DHomogeneous1D ( const LibUtilities::SessionReaderSharedPtr pSession,
const LibUtilities::BasisKey HomoBasis,
const NekDouble  lhom,
const bool  useFFT,
const bool  dealiasing,
const SpatialDomains::MeshGraphSharedPtr graph2D,
const std::string &  var = "DefaultVar",
const Collections::ImplementationType  ImpType = Collections::eNoImpType 
)

Sets up a list of local expansions based on an input mesh.

Definition at line 58 of file ExpList3DHomogeneous1D.cpp.

64 : ExpListHomogeneous1D(e3DH1D, pSession, HomoBasis, lhom, useFFT,
65 dealiasing)
66{
67 GenExpList3DHomogeneous1D(graph2D->GetExpansionInfo(var), ImpType);
69 m_graph = graph2D;
70}
void GenExpList3DHomogeneous1D(const SpatialDomains::ExpansionInfoMap &expansions, const Collections::ImplementationType ImpType)
SpatialDomains::MeshGraphSharedPtr m_graph
Mesh associated with this expansion list.
Definition: ExpList.h:1057
ExpansionType m_expType
Exapnsion type.
Definition: ExpList.h:1048

References Nektar::MultiRegions::e3DH1D, GenExpList3DHomogeneous1D(), Nektar::MultiRegions::ExpList::m_expType, and Nektar::MultiRegions::ExpList::m_graph.

◆ ExpList3DHomogeneous1D() [4/6]

Nektar::MultiRegions::ExpList3DHomogeneous1D::ExpList3DHomogeneous1D ( const LibUtilities::SessionReaderSharedPtr pSession,
const LibUtilities::BasisKey HomoBasis,
const NekDouble  lhom,
const bool  useFFT,
const bool  dealiasing,
const SpatialDomains::ExpansionInfoMap expansions,
const Collections::ImplementationType  ImpType = Collections::eNoImpType 
)

Sets up a list of local expansions based on an mesh expansion.

Definition at line 73 of file ExpList3DHomogeneous1D.cpp.

79 : ExpListHomogeneous1D(e3DH1D, pSession, HomoBasis, lhom, useFFT,
80 dealiasing)
81{
82 GenExpList3DHomogeneous1D(expansions, ImpType);
83}

References GenExpList3DHomogeneous1D().

◆ ExpList3DHomogeneous1D() [5/6]

Nektar::MultiRegions::ExpList3DHomogeneous1D::ExpList3DHomogeneous1D ( const ExpList3DHomogeneous1D In,
const bool  DeclarePlanesSetCoeffPhys = true 
)

Copy constructor.

Parameters
InExpList3DHomogeneous1D object to copy.

Definition at line 121 of file ExpList3DHomogeneous1D.cpp.

124{
125 if (DeclarePlanesSetCoeffPhys)
126 {
127 bool False = false;
128 ExpListSharedPtr zero_plane =
129 std::dynamic_pointer_cast<ExpList>(In.m_planes[0]);
130
131 for (int n = 0; n < m_planes.size(); ++n)
132 {
133 m_planes[n] =
135 }
136
137 SetCoeffPhys();
138 }
139}
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void SetCoeffPhys(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
Array< OneD, ExpListSharedPtr > m_planes
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::MultiRegions::ExpListHomogeneous1D::m_planes, and SetCoeffPhys().

◆ ExpList3DHomogeneous1D() [6/6]

Nektar::MultiRegions::ExpList3DHomogeneous1D::ExpList3DHomogeneous1D ( const ExpList3DHomogeneous1D In,
const std::vector< unsigned int > &  eIDs,
const bool  DeclarePlanesSetCoeffPhys = true,
const Collections::ImplementationType  ImpType = Collections::eNoImpType 
)

Constructor copying only elements defined in eIds.

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

Definition at line 145 of file ExpList3DHomogeneous1D.cpp.

149 : ExpListHomogeneous1D(In, eIDs, ImpType)
150{
151 if (DeclarePlanesSetCoeffPhys)
152 {
153 bool False = false;
154 std::vector<unsigned int> eIDsPlane;
155 int nel = eIDs.size() / m_planes.size();
156
157 for (int i = 0; i < nel; ++i)
158 {
159 eIDsPlane.push_back(eIDs[i]);
160 }
161
162 ExpListSharedPtr zero_plane_old =
163 std::dynamic_pointer_cast<ExpList>(In.m_planes[0]);
164
166 *(zero_plane_old), eIDsPlane, DeclarePlanesSetCoeffPhys, ImpType);
167
168 for (int n = 0; n < m_planes.size(); ++n)
169 {
170 m_planes[n] =
172 }
173
174 SetCoeffPhys();
175 }
176}

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::MultiRegions::ExpListHomogeneous1D::m_planes, and SetCoeffPhys().

◆ ~ExpList3DHomogeneous1D()

Nektar::MultiRegions::ExpList3DHomogeneous1D::~ExpList3DHomogeneous1D ( )
override

Destructor.

Destructor

Definition at line 181 of file ExpList3DHomogeneous1D.cpp.

182{
183}

Member Function Documentation

◆ GenExpList3DHomogeneous1D()

void Nektar::MultiRegions::ExpList3DHomogeneous1D::GenExpList3DHomogeneous1D ( const SpatialDomains::ExpansionInfoMap expansions,
const Collections::ImplementationType  ImpType 
)
private

Definition at line 85 of file ExpList3DHomogeneous1D.cpp.

88{
89 int n, j, nel;
90 bool False = false;
91 ExpListSharedPtr plane_zero;
92
93 // note that nzplanes can be larger than nzmodes
95 m_session, expansions, false, ImpType);
96
98 nel = m_planes[0]->GetExpSize();
99
100 for (j = 0; j < nel; ++j)
101 {
102 (*m_exp).push_back(m_planes[0]->GetExp(j));
103 }
104
105 for (n = 1; n < m_planes.size(); ++n)
106 {
107 m_planes[n] =
109 for (j = 0; j < nel; ++j)
110 {
111 (*m_exp).push_back((*m_exp)[j]);
112 }
113 }
114
115 SetCoeffPhys();
116}
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1115
const std::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:2079
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:1055

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::MultiRegions::ExpList::GetExp(), Nektar::MultiRegions::ExpList::m_exp, Nektar::MultiRegions::ExpListHomogeneous1D::m_planes, Nektar::MultiRegions::ExpList::m_session, and SetCoeffPhys().

Referenced by ExpList3DHomogeneous1D().

◆ SetCoeffPhys()

void Nektar::MultiRegions::ExpList3DHomogeneous1D::SetCoeffPhys ( void  )
protected

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

Definition at line 185 of file ExpList3DHomogeneous1D.cpp.

186{
187 int i, n, cnt;
188 int ncoeffs_per_plane = m_planes[0]->GetNcoeffs();
189 int npoints_per_plane = m_planes[0]->GetTotPoints();
190
191 int nzplanes = m_planes.size();
192
193 // Set total coefficients and points
194 m_ncoeffs = ncoeffs_per_plane * nzplanes;
195 m_npoints = npoints_per_plane * nzplanes;
196
197 m_coeffs = Array<OneD, NekDouble>(m_ncoeffs, 0.0);
198 m_phys = Array<OneD, NekDouble>(m_npoints, 0.0);
199
200 int nel = m_planes[0]->GetExpSize();
201 m_coeff_offset = Array<OneD, int>(nel * nzplanes);
202 m_phys_offset = Array<OneD, int>(nel * nzplanes);
203 Array<OneD, NekDouble> tmparray;
204
205 for (cnt = n = 0; n < nzplanes; ++n)
206 {
207 m_planes[n]->SetCoeffsArray(tmparray =
208 m_coeffs + ncoeffs_per_plane * n);
209 m_planes[n]->SetPhysArray(tmparray = m_phys + npoints_per_plane * n);
210
211 for (i = 0; i < nel; ++i)
212 {
213 m_coeff_offset[cnt] =
214 m_planes[n]->GetCoeff_Offset(i) + n * ncoeffs_per_plane;
215 m_phys_offset[cnt++] =
216 m_planes[n]->GetPhys_Offset(i) + n * npoints_per_plane;
217 }
218 }
219}
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:1080
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:1124
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:1060
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:1126
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:1096

References Nektar::MultiRegions::ExpList::m_coeff_offset, Nektar::MultiRegions::ExpList::m_coeffs, Nektar::MultiRegions::ExpList::m_ncoeffs, Nektar::MultiRegions::ExpList::m_npoints, Nektar::MultiRegions::ExpList::m_phys, Nektar::MultiRegions::ExpList::m_phys_offset, and Nektar::MultiRegions::ExpListHomogeneous1D::m_planes.

Referenced by Nektar::MultiRegions::ContField3DHomogeneous1D::ContField3DHomogeneous1D(), Nektar::MultiRegions::DisContField3DHomogeneous1D::DisContField3DHomogeneous1D(), ExpList3DHomogeneous1D(), and GenExpList3DHomogeneous1D().

◆ v_GetCoords() [1/2]

void Nektar::MultiRegions::ExpList3DHomogeneous1D::v_GetCoords ( Array< OneD, NekDouble > &  xc0,
Array< OneD, NekDouble > &  xc1,
Array< OneD, NekDouble > &  xc2 
)
overrideprotectedvirtual

The operation calls the 2D plane coordinates through the function ExpList::GetCoords and then evaluated the third coordinate using the member m_lhom

Parameters
coord_0After calculation, the \(x_1\) coordinate will be stored in this array.
coord_1After calculation, the \(x_2\) coordinate will be stored in this array.
coord_2After calculation, the \(x_3\) coordinate will be stored in this array. This coordinate is evaluated using the predefined value m_lhom

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 275 of file ExpList3DHomogeneous1D.cpp.

278{
279 int n;
280 Array<OneD, NekDouble> tmp_xc;
281 int nzplanes = m_planes.size();
282 int npoints = m_planes[0]->GetTotPoints();
283
284 m_planes[0]->GetCoords(xc0, xc1);
285
286 // Fill z-direction
287 Array<OneD, const NekDouble> pts = m_homogeneousBasis->GetZ();
288
289 Array<OneD, NekDouble> local_pts(m_planes.size());
290
291 for (n = 0; n < m_planes.size(); n++)
292 {
293 local_pts[n] = pts[m_transposition->GetPlaneID(n)];
294 }
295
296 Array<OneD, NekDouble> z(nzplanes);
297
298 Vmath::Smul(nzplanes, m_lhom / 2.0, local_pts, 1, z, 1);
299 Vmath::Sadd(nzplanes, m_lhom / 2.0, z, 1, z, 1);
300
301 for (n = 0; n < nzplanes; ++n)
302 {
303 Vmath::Fill(npoints, z[n], tmp_xc = xc2 + npoints * n, 1);
304 if (n)
305 {
306 Vmath::Vcopy(npoints, xc0, 1, tmp_xc = xc0 + npoints * n, 1);
307 Vmath::Vcopy(npoints, xc1, 1, tmp_xc = xc1 + npoints * n, 1);
308 }
309 }
310}
NekDouble m_lhom
Width of homogeneous direction.
LibUtilities::TranspositionSharedPtr m_transposition
LibUtilities::BasisSharedPtr m_homogeneousBasis
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
std::vector< double > z(NPUPPER)
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.hpp:54
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition: Vmath.hpp:194
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References Vmath::Fill(), Nektar::MultiRegions::ExpListHomogeneous1D::m_homogeneousBasis, Nektar::MultiRegions::ExpListHomogeneous1D::m_lhom, Nektar::MultiRegions::ExpListHomogeneous1D::m_planes, Nektar::MultiRegions::ExpListHomogeneous1D::m_transposition, Vmath::Sadd(), Vmath::Smul(), Vmath::Vcopy(), and Nektar::UnitTests::z().

Referenced by v_WriteVtkPieceHeader().

◆ v_GetCoords() [2/2]

void Nektar::MultiRegions::ExpList3DHomogeneous1D::v_GetCoords ( const int  eid,
Array< OneD, NekDouble > &  xc0,
Array< OneD, NekDouble > &  xc1,
Array< OneD, NekDouble > &  xc2 
)
finalprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 222 of file ExpList3DHomogeneous1D.cpp.

226{
227 int n;
228 Array<OneD, NekDouble> tmp_xc;
229 int nzplanes = m_planes.size();
230 int npoints = GetTotPoints(eid);
231
232 (*m_exp)[eid]->GetCoords(xc0, xc1);
233
234 // Fill z-direction
235 Array<OneD, const NekDouble> pts = m_homogeneousBasis->GetZ();
236 Array<OneD, NekDouble> local_pts(m_planes.size());
237
238 for (n = 0; n < m_planes.size(); n++)
239 {
240 local_pts[n] = pts[m_transposition->GetPlaneID(n)];
241 }
242
243 Array<OneD, NekDouble> z(nzplanes);
244
245 Vmath::Smul(nzplanes, m_lhom / 2.0, local_pts, 1, z, 1);
246 Vmath::Sadd(nzplanes, m_lhom / 2.0, z, 1, z, 1);
247
248 for (n = 0; n < nzplanes; ++n)
249 {
250 Vmath::Fill(npoints, z[n], tmp_xc = xc2 + npoints * n, 1);
251 if (n)
252 {
253 Vmath::Vcopy(npoints, xc0, 1, tmp_xc = xc0 + npoints * n, 1);
254 Vmath::Vcopy(npoints, xc1, 1, tmp_xc = xc1 + npoints * n, 1);
255 }
256 }
257}
int GetTotPoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1557

References Vmath::Fill(), Nektar::MultiRegions::ExpList::GetTotPoints(), Nektar::MultiRegions::ExpListHomogeneous1D::m_homogeneousBasis, Nektar::MultiRegions::ExpListHomogeneous1D::m_lhom, Nektar::MultiRegions::ExpListHomogeneous1D::m_planes, Nektar::MultiRegions::ExpListHomogeneous1D::m_transposition, Vmath::Sadd(), Vmath::Smul(), Vmath::Vcopy(), and Nektar::UnitTests::z().

◆ v_GetPeriodicEntities()

void Nektar::MultiRegions::ExpList3DHomogeneous1D::v_GetPeriodicEntities ( PeriodicMap periodicVerts,
PeriodicMap periodicEdges,
PeriodicMap periodicFaces 
)
inlineoverrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 129 of file ExpList3DHomogeneous1D.h.

132 {
133 m_planes[0]->GetPeriodicEntities(periodicVerts, periodicEdges);
134 }

References Nektar::MultiRegions::ExpListHomogeneous1D::m_planes.

◆ v_HomogeneousEnergy()

Array< OneD, const NekDouble > Nektar::MultiRegions::ExpList3DHomogeneous1D::v_HomogeneousEnergy ( void  )
overrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 501 of file ExpList3DHomogeneous1D.cpp.

502{
503 Array<OneD, NekDouble> energy(m_planes.size() / 2);
504 NekDouble area = 0.0;
505
506 // Calculate total area of elements.
507 for (int n = 0; n < m_planes[0]->GetExpSize(); ++n)
508 {
509 Array<OneD, NekDouble> inarray(m_planes[0]->GetExp(n)->GetTotPoints(),
510 1.0);
511 area += m_planes[0]->GetExp(n)->Integral(inarray);
512 }
513
514 m_comm->GetRowComm()->AllReduce(area, LibUtilities::ReduceSum);
515
516 // Calculate L2 norm of real/imaginary planes.
517 for (int n = 0; n < m_planes.size(); n += 2)
518 {
519 NekDouble err;
520
521 energy[n / 2] = 0;
522
523 for (int i = 0; i < m_planes[n]->GetExpSize(); ++i)
524 {
525 LocalRegions::ExpansionSharedPtr exp = m_planes[n]->GetExp(i);
526 Array<OneD, NekDouble> phys(exp->GetTotPoints());
527 exp->BwdTrans(m_planes[n]->GetCoeffs() +
529 phys);
530 err = exp->L2(phys);
531 energy[n / 2] += err * err;
532
533 exp = m_planes[n + 1]->GetExp(i);
534 exp->BwdTrans(m_planes[n + 1]->GetCoeffs() +
535 m_planes[n + 1]->GetCoeff_Offset(i),
536 phys);
537 err = exp->L2(phys);
538 energy[n / 2] += err * err;
539 }
540
541 m_comm->GetRowComm()->AllReduce(energy[n / 2], LibUtilities::ReduceSum);
542 energy[n / 2] /= 2.0 * area;
543 }
544
545 return energy;
546}
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:1952
int GetCoeff_Offset(int n) const
Get the start offset position for a local contiguous list of coeffs correspoinding to element n.
Definition: ExpList.h:2087
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: ExpList.h:1053
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:66
double NekDouble

References Nektar::MultiRegions::ExpList::GetCoeff_Offset(), Nektar::MultiRegions::ExpList::GetCoeffs(), Nektar::MultiRegions::ExpList::GetExp(), Nektar::MultiRegions::ExpList::GetTotPoints(), Nektar::MultiRegions::ExpList::m_comm, Nektar::MultiRegions::ExpListHomogeneous1D::m_planes, and Nektar::LibUtilities::ReduceSum.

◆ v_L2()

NekDouble Nektar::MultiRegions::ExpList3DHomogeneous1D::v_L2 ( const Array< OneD, const NekDouble > &  inarray,
const Array< OneD, const NekDouble > &  soln = NullNekDouble1DArray 
)
overrideprotectedvirtual

Given a spectral/hp approximation \(u^{\delta}(\boldsymbol{x})\) evaluated at the quadrature points (which should be contained in m_phys), this function calculates the \(L_2\) error of this approximation with respect to an exact solution. The local distribution of the quadrature points allows an elemental evaluation of this operation through the functions StdRegions::StdExpansion::L2.

The exact solution, also evaluated at the quadrature points, should be contained in the variable m_phys of the ExpList object Sol.

Parameters
SolAn ExpList, containing the discrete evaluation of the exact solution at the quadrature points in its array m_phys.
Returns
The \(L_2\) error of the approximation.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 464 of file ExpList3DHomogeneous1D.cpp.

467{
468 int cnt = 0;
469 NekDouble errL2, err = 0.0;
470 Array<OneD, const NekDouble> w = m_homogeneousBasis->GetW();
471 Array<OneD, NekDouble> local_w(m_planes.size());
472
473 for (int n = 0; n < m_planes.size(); n++)
474 {
475 local_w[n] = w[m_transposition->GetPlaneID(n)];
476 }
477
478 if (soln == NullNekDouble1DArray)
479 {
480 for (int n = 0; n < m_planes.size(); ++n)
481 {
482 errL2 = m_planes[n]->L2(inarray + cnt);
483 cnt += m_planes[n]->GetTotPoints();
484 err += errL2 * errL2 * local_w[n] * m_lhom * 0.5;
485 }
486 }
487 else
488 {
489 for (int n = 0; n < m_planes.size(); ++n)
490 {
491 errL2 = m_planes[n]->L2(inarray + cnt, soln + cnt);
492 cnt += m_planes[n]->GetTotPoints();
493 err += errL2 * errL2 * local_w[n] * m_lhom * 0.5;
494 }
495 }
496 m_comm->GetColumnComm()->AllReduce(err, LibUtilities::ReduceSum);
497
498 return sqrt(err);
499}
std::vector< double > w(NPUPPER)
static Array< OneD, NekDouble > NullNekDouble1DArray
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References Nektar::MultiRegions::ExpList::m_comm, Nektar::MultiRegions::ExpListHomogeneous1D::m_homogeneousBasis, Nektar::MultiRegions::ExpListHomogeneous1D::m_lhom, Nektar::MultiRegions::ExpListHomogeneous1D::m_planes, Nektar::MultiRegions::ExpListHomogeneous1D::m_transposition, Nektar::NullNekDouble1DArray, Nektar::LibUtilities::ReduceSum, tinysimd::sqrt(), and Nektar::UnitTests::w().

◆ v_WriteTecplotConnectivity()

void Nektar::MultiRegions::ExpList3DHomogeneous1D::v_WriteTecplotConnectivity ( std::ostream &  outfile,
int  expansion 
)
overrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 312 of file ExpList3DHomogeneous1D.cpp.

314{
315 ASSERTL0(expansion == -1,
316 "Multi-zone output not supported for homogeneous expansions.");
317
318 const int nPtsPlane = m_planes[0]->GetNpoints();
319 const int nElmt = m_planes[0]->GetExpSize();
320 const int nPlanes = m_planes.size();
321
322 for (int i = 0, cnt = 0; i < nElmt; ++i)
323 {
324 const int np0 = (*m_exp)[i]->GetNumPoints(0);
325 const int np1 = (*m_exp)[i]->GetNumPoints(1);
326
327 for (int n = 1; n < nPlanes; ++n)
328 {
329 const int o1 = (n - 1) * nPtsPlane;
330 const int o2 = n * nPtsPlane;
331 for (int j = 1; j < np1; ++j)
332 {
333 for (int k = 1; k < np0; ++k)
334 {
335 outfile << cnt + (j - 1) * np0 + (k - 1) + o1 + 1 << " ";
336 outfile << cnt + (j - 1) * np0 + (k - 1) + o2 + 1 << " ";
337 outfile << cnt + (j - 1) * np0 + k + o2 + 1 << " ";
338 outfile << cnt + (j - 1) * np0 + k + o1 + 1 << " ";
339 outfile << cnt + j * np0 + (k - 1) + o1 + 1 << " ";
340 outfile << cnt + j * np0 + (k - 1) + o2 + 1 << " ";
341 outfile << cnt + j * np0 + k + o2 + 1 << " ";
342 outfile << cnt + j * np0 + k + o1 + 1 << endl;
343 }
344 }
345 }
346
347 cnt += np0 * np1;
348 }
349}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208

References ASSERTL0, and Nektar::MultiRegions::ExpListHomogeneous1D::m_planes.

◆ v_WriteVtkPieceHeader()

void Nektar::MultiRegions::ExpList3DHomogeneous1D::v_WriteVtkPieceHeader ( std::ostream &  outfile,
int  expansion,
int  istrip 
)
overrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 351 of file ExpList3DHomogeneous1D.cpp.

353{
354 // If there is only one plane (e.g. HalfMode), we write a 2D plane.
355 if (m_planes.size() == 1)
356 {
357 m_planes[0]->WriteVtkPieceHeader(outfile, expansion);
358 return;
359 }
360
361 // If we are using Fourier points, output extra plane to fill domain
362 int outputExtraPlane = 0;
363 if (m_homogeneousBasis->GetBasisType() == LibUtilities::eFourier &&
364 m_homogeneousBasis->GetPointsType() ==
366 {
367 outputExtraPlane = 1;
368 }
369
370 int i, j, k;
371 int nq0 = (*m_exp)[expansion]->GetNumPoints(0);
372 int nq1 = (*m_exp)[expansion]->GetNumPoints(1);
373 int nq2 = m_planes.size() + outputExtraPlane;
374 int ntot = nq0 * nq1 * nq2;
375 int ntotminus = (nq0 - 1) * (nq1 - 1) * (nq2 - 1);
376
377 Array<OneD, NekDouble> coords[3];
378 coords[0] = Array<OneD, NekDouble>(ntot);
379 coords[1] = Array<OneD, NekDouble>(ntot);
380 coords[2] = Array<OneD, NekDouble>(ntot);
381 v_GetCoords(expansion, coords[0], coords[1], coords[2]);
382
383 if (outputExtraPlane)
384 {
385 // Copy coords[0] and coords[1] to extra plane
386 Array<OneD, NekDouble> tmp;
387 Vmath::Vcopy(nq0 * nq1, coords[0], 1,
388 tmp = coords[0] + (nq2 - 1) * nq0 * nq1, 1);
389 Vmath::Vcopy(nq0 * nq1, coords[1], 1,
390 tmp = coords[1] + (nq2 - 1) * nq0 * nq1, 1);
391 // Fill coords[2] for extra plane
392 NekDouble z = coords[2][nq0 * nq1 * m_planes.size() - 1] +
393 (coords[2][nq0 * nq1] - coords[2][0]);
394 Vmath::Fill(nq0 * nq1, z, tmp = coords[2] + (nq2 - 1) * nq0 * nq1, 1);
395 }
396
397 NekDouble DistStrip;
398 m_session->LoadParameter("DistStrip", DistStrip, 0);
399 // Reset the z-coords for homostrips
400 for (int i = 0; i < ntot; i++)
401 {
402 coords[2][i] += istrip * DistStrip;
403 }
404
405 outfile << " <Piece NumberOfPoints=\"" << ntot << "\" NumberOfCells=\""
406 << ntotminus << "\">" << endl;
407 outfile << " <Points>" << endl;
408 outfile << " <DataArray type=\"Float64\" "
409 << R"(NumberOfComponents="3" format="ascii">)" << endl;
410 outfile << " ";
411 for (i = 0; i < ntot; ++i)
412 {
413 for (j = 0; j < 3; ++j)
414 {
415 outfile << coords[j][i] << " ";
416 }
417 outfile << endl;
418 }
419 outfile << endl;
420 outfile << " </DataArray>" << endl;
421 outfile << " </Points>" << endl;
422 outfile << " <Cells>" << endl;
423 outfile << " <DataArray type=\"Int32\" "
424 << R"(Name="connectivity" format="ascii">)" << endl;
425 for (i = 0; i < nq0 - 1; ++i)
426 {
427 for (j = 0; j < nq1 - 1; ++j)
428 {
429 for (k = 0; k < nq2 - 1; ++k)
430 {
431 outfile << k * nq0 * nq1 + j * nq0 + i << " "
432 << k * nq0 * nq1 + j * nq0 + i + 1 << " "
433 << k * nq0 * nq1 + (j + 1) * nq0 + i + 1 << " "
434 << k * nq0 * nq1 + (j + 1) * nq0 + i << " "
435 << (k + 1) * nq0 * nq1 + j * nq0 + i << " "
436 << (k + 1) * nq0 * nq1 + j * nq0 + i + 1 << " "
437 << (k + 1) * nq0 * nq1 + (j + 1) * nq0 + i + 1 << " "
438 << (k + 1) * nq0 * nq1 + (j + 1) * nq0 + i << endl;
439 }
440 }
441 }
442 outfile << endl;
443 outfile << " </DataArray>" << endl;
444 outfile << " <DataArray type=\"Int32\" "
445 << R"(Name="offsets" format="ascii">)" << endl;
446 for (i = 0; i < ntotminus; ++i)
447 {
448 outfile << i * 8 + 8 << " ";
449 }
450 outfile << endl;
451 outfile << " </DataArray>" << endl;
452 outfile << " <DataArray type=\"UInt8\" "
453 << R"(Name="types" format="ascii">)" << endl;
454 for (i = 0; i < ntotminus; ++i)
455 {
456 outfile << "12 ";
457 }
458 outfile << endl;
459 outfile << " </DataArray>" << endl;
460 outfile << " </Cells>" << endl;
461 outfile << " <PointData>" << endl;
462}
void v_GetCoords(Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1, Array< OneD, NekDouble > &coord_2) override
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:74
@ eFourier
Fourier Expansion .
Definition: BasisType.h:55

References Nektar::LibUtilities::eFourier, Nektar::LibUtilities::eFourierEvenlySpaced, Vmath::Fill(), Nektar::MultiRegions::ExpListHomogeneous1D::m_homogeneousBasis, Nektar::MultiRegions::ExpListHomogeneous1D::m_planes, Nektar::MultiRegions::ExpList::m_session, v_GetCoords(), Vmath::Vcopy(), and Nektar::UnitTests::z().