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

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

#include <ExpListHomogeneous2D.h>

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

Public Member Functions

 ExpListHomogeneous2D (const ExpansionType type)
 Default constructor. More...
 
 ExpListHomogeneous2D (const ExpansionType type, const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::BasisKey &HomoBasis_y, const LibUtilities::BasisKey &HomoBasis_z, const NekDouble ly, const NekDouble lz, const bool useFFT, const bool dealiasing)
 
 ExpListHomogeneous2D (const ExpListHomogeneous2D &In)
 Copy constructor. More...
 
 ExpListHomogeneous2D (const ExpListHomogeneous2D &In, const std::vector< unsigned int > &eIDs)
 
 ~ExpListHomogeneous2D () override
 Destructor. More...
 
void Homogeneous2DTrans (const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool IsForwards, bool Shuff=true, bool UnShuff=true)
 
void SetPaddingBase (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)
 
- 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...
 

Public Attributes

bool m_useFFT
 FFT variables. More...
 
LibUtilities::NektarFFTSharedPtr m_FFT_y
 
LibUtilities::NektarFFTSharedPtr m_FFT_z
 
Array< OneD, NekDoublem_tmpIN
 
Array< OneD, NekDoublem_tmpOUT
 
LibUtilities::TranspositionSharedPtr m_transposition
 
LibUtilities::CommSharedPtr m_Ycomm
 
LibUtilities::CommSharedPtr m_Zcomm
 

Protected Member Functions

DNekBlkMatSharedPtr GenHomogeneous2DBlockMatrix (Homogeneous2DMatType mattype) const
 
DNekBlkMatSharedPtr GetHomogeneous2DBlockMatrix (Homogeneous2DMatType mattype) const
 
size_t v_GetNumElmts (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
 
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_WriteVtkPieceData (std::ostream &outfile, int expansion, std::string var) 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
 
- 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)
 

Protected Attributes

LibUtilities::BasisSharedPtr m_homogeneousBasis_y
 Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m_coeff and m_phys. More...
 
LibUtilities::BasisSharedPtr m_homogeneousBasis_z
 Base expansion in z direction. More...
 
LibUtilities::BasisSharedPtr m_paddingBasis_y
 Base expansion in y direction. More...
 
LibUtilities::BasisSharedPtr m_paddingBasis_z
 Base expansion in z direction. More...
 
NekDouble m_lhom_y
 Width of homogeneous direction y. More...
 
NekDouble m_lhom_z
 Width of homogeneous direction z. More...
 
Homo2DBlockMatrixMapShPtr m_homogeneous2DBlockMat
 
Array< OneD, ExpListSharedPtrm_lines
 Vector of ExpList, will be filled with ExpList1D. More...
 
int m_ny
 Number of modes = number of poitns in y direction. More...
 
int m_nz
 Number of modes = number of poitns in z direction. More...
 
- 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...
 

Private Attributes

bool m_dealiasing
 
int m_padsize_y
 
int m_padsize_z
 
DNekMatSharedPtr MatFwdPAD
 
DNekMatSharedPtr MatBwdPAD
 

Additional Inherited Members

- Static Protected Member Functions inherited from Nektar::MultiRegions::ExpList
static SpatialDomains::BoundaryConditionShPtr GetBoundaryCondition (const SpatialDomains::BoundaryConditionCollection &collection, unsigned int index, const std::string &variable)
 

Detailed Description

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

Definition at line 77 of file ExpListHomogeneous2D.h.

Constructor & Destructor Documentation

◆ ExpListHomogeneous2D() [1/4]

Nektar::MultiRegions::ExpListHomogeneous2D::ExpListHomogeneous2D ( const ExpansionType  type)

Default constructor.

Definition at line 46 of file ExpListHomogeneous2D.cpp.

51{
52}
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
LibUtilities::BasisSharedPtr m_homogeneousBasis_y
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
NekDouble m_lhom_z
Width of homogeneous direction z.
LibUtilities::BasisSharedPtr m_homogeneousBasis_z
Base expansion in z direction.
NekDouble m_lhom_y
Width of homogeneous direction y.
ExpList(const ExpansionType Type=eNoType)
The default constructor using a type.
Definition: ExpList.cpp:90
static BasisSharedPtr NullBasisSharedPtr
Definition: Basis.h:350

◆ ExpListHomogeneous2D() [2/4]

Nektar::MultiRegions::ExpListHomogeneous2D::ExpListHomogeneous2D ( const ExpansionType  type,
const LibUtilities::SessionReaderSharedPtr pSession,
const LibUtilities::BasisKey HomoBasis_y,
const LibUtilities::BasisKey HomoBasis_z,
const NekDouble  ly,
const NekDouble  lz,
const bool  useFFT,
const bool  dealiasing 
)

Definition at line 54 of file ExpListHomogeneous2D.cpp.

60 : ExpList(type), m_useFFT(useFFT), m_lhom_y(lhom_y), m_lhom_z(lhom_z),
63 m_dealiasing(dealiasing)
64{
65 m_session = pSession;
66 m_comm = pSession->GetComm();
67
69 "Homogeneous Basis in y direction is a null basis");
71 "Homogeneous Basis in z direction is a null basis");
72
75
78 HomoBasis_y, HomoBasis_z, m_comm->GetColumnComm());
79
80 m_Ycomm = m_comm->GetColumnComm()->GetRowComm();
81 m_Zcomm = m_comm->GetColumnComm()->GetRowComm();
82
83 m_ny = m_homogeneousBasis_y->GetNumPoints() / m_Ycomm->GetSize();
84 m_nz = m_homogeneousBasis_z->GetNumPoints() / m_Zcomm->GetSize();
85
86 m_lines = Array<OneD, ExpListSharedPtr>(m_ny * m_nz);
87
88 if (m_useFFT)
89 {
90 m_FFT_y =
92 m_FFT_z =
94 }
95
96 if (m_dealiasing)
97 {
98 ASSERTL0(m_comm->GetColumnComm()->GetSize() == 1,
99 "Remove dealiasing if you want to run in parallel");
101 }
102}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:265
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:143
LibUtilities::NektarFFTSharedPtr m_FFT_y
int m_nz
Number of modes = number of poitns in z direction.
LibUtilities::NektarFFTSharedPtr m_FFT_z
int m_ny
Number of modes = number of poitns in y direction.
LibUtilities::TranspositionSharedPtr m_transposition
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: ExpList.h:1053
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:1055
BasisManagerT & BasisManager(void)
NektarFFTFactory & GetNektarFFTFactory()
Definition: NektarFFT.cpp:65
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, ASSERTL2, Nektar::LibUtilities::BasisManager(), Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::GetNektarFFTFactory(), Nektar::MultiRegions::ExpList::m_comm, m_dealiasing, m_FFT_y, m_FFT_z, m_homogeneousBasis_y, m_homogeneousBasis_z, m_lines, m_ny, m_nz, Nektar::MultiRegions::ExpList::m_session, m_transposition, m_useFFT, m_Ycomm, m_Zcomm, Nektar::LibUtilities::NullBasisKey(), and SetPaddingBase().

◆ ExpListHomogeneous2D() [3/4]

Nektar::MultiRegions::ExpListHomogeneous2D::ExpListHomogeneous2D ( const ExpListHomogeneous2D In)

Copy constructor.

Parameters
InExpListHomogeneous2D object to copy.

Definition at line 107 of file ExpListHomogeneous2D.cpp.

108 : ExpList(In, false), m_useFFT(In.m_useFFT), m_FFT_y(In.m_FFT_y),
109 m_FFT_z(In.m_FFT_z), m_transposition(In.m_transposition),
110 m_Ycomm(In.m_Ycomm), m_Zcomm(In.m_Ycomm),
111 m_homogeneousBasis_y(In.m_homogeneousBasis_y),
112 m_homogeneousBasis_z(In.m_homogeneousBasis_z), m_lhom_y(In.m_lhom_y),
113 m_lhom_z(In.m_lhom_z),
114 m_homogeneous2DBlockMat(In.m_homogeneous2DBlockMat), m_ny(In.m_ny),
115 m_nz(In.m_nz), m_dealiasing(In.m_dealiasing), m_padsize_y(In.m_padsize_y),
116 m_padsize_z(In.m_padsize_z), MatFwdPAD(In.MatFwdPAD),
117 MatBwdPAD(In.MatBwdPAD)
118{
119 m_lines = Array<OneD, ExpListSharedPtr>(In.m_lines.size());
120}

References m_lines.

◆ ExpListHomogeneous2D() [4/4]

Nektar::MultiRegions::ExpListHomogeneous2D::ExpListHomogeneous2D ( const ExpListHomogeneous2D In,
const std::vector< unsigned int > &  eIDs 
)

Definition at line 122 of file ExpListHomogeneous2D.cpp.

124 : ExpList(In, eIDs, false), m_useFFT(In.m_useFFT), m_FFT_y(In.m_FFT_y),
125 m_FFT_z(In.m_FFT_z), m_transposition(In.m_transposition),
126 m_Ycomm(In.m_Ycomm), m_Zcomm(In.m_Ycomm),
127 m_homogeneousBasis_y(In.m_homogeneousBasis_y),
128 m_homogeneousBasis_z(In.m_homogeneousBasis_z), m_lhom_y(In.m_lhom_y),
129 m_lhom_z(In.m_lhom_z),
132 m_ny(In.m_ny), m_nz(In.m_nz), m_dealiasing(In.m_dealiasing),
133 m_padsize_y(In.m_padsize_y), m_padsize_z(In.m_padsize_z),
134 MatFwdPAD(In.MatFwdPAD), MatBwdPAD(In.MatBwdPAD)
135{
136 m_lines = Array<OneD, ExpListSharedPtr>(In.m_lines.size());
137}

References m_lines.

◆ ~ExpListHomogeneous2D()

Nektar::MultiRegions::ExpListHomogeneous2D::~ExpListHomogeneous2D ( )
override

Destructor.

Destructor

Definition at line 142 of file ExpListHomogeneous2D.cpp.

143{
144}

Member Function Documentation

◆ GenHomogeneous2DBlockMatrix()

DNekBlkMatSharedPtr Nektar::MultiRegions::ExpListHomogeneous2D::GenHomogeneous2DBlockMatrix ( Homogeneous2DMatType  mattype) const
protected

Definition at line 544 of file ExpListHomogeneous2D.cpp.

546{
547 int i;
548 int n_exp = 0;
549
550 DNekMatSharedPtr loc_mat;
551 DNekBlkMatSharedPtr BlkMatrix;
552
554
555 int NumPoints = 0;
556 int NumModes = 0;
557 int NumPencils = 0;
558
559 if ((mattype == eForwardsCoeffSpaceY1D) ||
560 (mattype == eBackwardsCoeffSpaceY1D) ||
561 (mattype == eForwardsPhysSpaceY1D) ||
562 (mattype == eBackwardsPhysSpaceY1D))
563 {
564 Basis = m_homogeneousBasis_y;
565 NumPoints = m_homogeneousBasis_y->GetNumModes();
566 NumModes = m_homogeneousBasis_y->GetNumPoints();
567 NumPencils = m_homogeneousBasis_z->GetNumPoints();
568 }
569 else
570 {
571 Basis = m_homogeneousBasis_z;
572 NumPoints = m_homogeneousBasis_z->GetNumModes();
573 NumModes = m_homogeneousBasis_z->GetNumPoints();
574 NumPencils = m_homogeneousBasis_y->GetNumPoints();
575 }
576
577 if ((mattype == eForwardsCoeffSpaceY1D) ||
578 (mattype == eForwardsCoeffSpaceZ1D) ||
579 (mattype == eBackwardsCoeffSpaceY1D) ||
580 (mattype == eBackwardsCoeffSpaceZ1D))
581 {
582 n_exp = m_lines[0]->GetNcoeffs();
583 }
584 else
585 {
586 n_exp = m_lines[0]->GetTotPoints(); // will operatore on m_phys
587 }
588
589 Array<OneD, unsigned int> nrows(n_exp);
590 Array<OneD, unsigned int> ncols(n_exp);
591
592 if ((mattype == eForwardsCoeffSpaceY1D) ||
593 (mattype == eForwardsPhysSpaceY1D) ||
594 (mattype == eForwardsCoeffSpaceZ1D) ||
595 (mattype == eForwardsPhysSpaceZ1D))
596 {
597 nrows = Array<OneD, unsigned int>(n_exp * NumPencils, NumModes);
598 ncols = Array<OneD, unsigned int>(n_exp * NumPencils, NumPoints);
599 }
600 else
601 {
602 nrows = Array<OneD, unsigned int>(n_exp * NumPencils, NumPoints);
603 ncols = Array<OneD, unsigned int>(n_exp * NumPencils, NumModes);
604 }
605
606 MatrixStorage blkmatStorage = eDIAGONAL;
607 BlkMatrix = MemoryManager<DNekBlkMat>::AllocateSharedPtr(nrows, ncols,
608 blkmatStorage);
609
610 StdRegions::StdSegExp StdSeg(Basis->GetBasisKey());
611
612 if ((mattype == eForwardsCoeffSpaceY1D) ||
613 (mattype == eForwardsPhysSpaceY1D) ||
614 (mattype == eForwardsCoeffSpaceZ1D) ||
615 (mattype == eForwardsPhysSpaceZ1D))
616 {
617 StdRegions::StdMatrixKey matkey(StdRegions::eFwdTrans,
618 StdSeg.DetShapeType(), StdSeg);
619
620 loc_mat = StdSeg.GetStdMatrix(matkey);
621 }
622 else
623 {
624 StdRegions::StdMatrixKey matkey(StdRegions::eBwdTrans,
625 StdSeg.DetShapeType(), StdSeg);
626
627 loc_mat = StdSeg.GetStdMatrix(matkey);
628 }
629
630 // set up array of block matrices.
631 for (i = 0; i < (n_exp * NumPencils); ++i)
632 {
633 BlkMatrix->SetBlock(i, i, loc_mat);
634 }
635
636 return BlkMatrix;
637}
std::shared_ptr< Basis > BasisSharedPtr
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:77
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::MultiRegions::eBackwardsCoeffSpaceY1D, Nektar::MultiRegions::eBackwardsCoeffSpaceZ1D, Nektar::MultiRegions::eBackwardsPhysSpaceY1D, Nektar::StdRegions::eBwdTrans, Nektar::eDIAGONAL, Nektar::MultiRegions::eForwardsCoeffSpaceY1D, Nektar::MultiRegions::eForwardsCoeffSpaceZ1D, Nektar::MultiRegions::eForwardsPhysSpaceY1D, Nektar::MultiRegions::eForwardsPhysSpaceZ1D, Nektar::StdRegions::eFwdTrans, Nektar::LibUtilities::Basis::GetBasisKey(), Nektar::StdRegions::StdExpansion::GetStdMatrix(), m_homogeneousBasis_y, m_homogeneousBasis_z, and m_lines.

Referenced by GetHomogeneous2DBlockMatrix().

◆ GetHomogeneous2DBlockMatrix()

DNekBlkMatSharedPtr Nektar::MultiRegions::ExpListHomogeneous2D::GetHomogeneous2DBlockMatrix ( Homogeneous2DMatType  mattype) const
protected

Definition at line 528 of file ExpListHomogeneous2D.cpp.

530{
531 auto matrixIter = m_homogeneous2DBlockMat->find(mattype);
532
533 if (matrixIter == m_homogeneous2DBlockMat->end())
534 {
535 return ((*m_homogeneous2DBlockMat)[mattype] =
537 }
538 else
539 {
540 return matrixIter->second;
541 }
542}
DNekBlkMatSharedPtr GenHomogeneous2DBlockMatrix(Homogeneous2DMatType mattype) const

References GenHomogeneous2DBlockMatrix(), and m_homogeneous2DBlockMat.

Referenced by Homogeneous2DTrans().

◆ Homogeneous2DTrans()

void Nektar::MultiRegions::ExpListHomogeneous2D::Homogeneous2DTrans ( const int  npts,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
bool  IsForwards,
bool  Shuff = true,
bool  UnShuff = true 
)

Definition at line 394 of file ExpListHomogeneous2D.cpp.

398{
399 if (m_useFFT)
400 {
401
402 int n = m_lines.size(); // number of Fourier points in the Fourier
403 // directions (x-z grid)
404 int p =
405 num_dofs /
406 n; // number of points per line = n of Fourier transform required
407
408 Array<OneD, NekDouble> fft_in(num_dofs);
409 Array<OneD, NekDouble> fft_out(num_dofs);
410
411 m_transposition->Transpose(num_dofs, inarray, fft_in, false,
413
414 if (IsForwards)
415 {
416 for (int i = 0; i < (p * m_nz); i++)
417 {
418 m_FFT_y->FFTFwdTrans(m_tmpIN = fft_in + i * m_ny,
419 m_tmpOUT = fft_out + i * m_ny);
420 }
421 }
422 else
423 {
424 for (int i = 0; i < (p * m_nz); i++)
425 {
426 m_FFT_y->FFTBwdTrans(m_tmpIN = fft_in + i * m_ny,
427 m_tmpOUT = fft_out + i * m_ny);
428 }
429 }
430
431 m_transposition->Transpose(num_dofs, fft_out, fft_in, false,
433
434 if (IsForwards)
435 {
436 for (int i = 0; i < (p * m_ny); i++)
437 {
438 m_FFT_z->FFTFwdTrans(m_tmpIN = fft_in + i * m_nz,
439 m_tmpOUT = fft_out + i * m_nz);
440 }
441 }
442 else
443 {
444 for (int i = 0; i < (p * m_ny); i++)
445 {
446 m_FFT_z->FFTBwdTrans(m_tmpIN = fft_in + i * m_nz,
447 m_tmpOUT = fft_out + i * m_nz);
448 }
449 }
450
451 // TODO: required ZYtoX routine
452 m_transposition->Transpose(num_dofs, fft_out, fft_in, false,
454
455 m_transposition->Transpose(num_dofs, fft_in, outarray, false,
457 }
458 else
459 {
460 DNekBlkMatSharedPtr blkmatY;
461 DNekBlkMatSharedPtr blkmatZ;
462
463 if (num_dofs == m_npoints) // transform phys space
464 {
465 if (IsForwards)
466 {
469 }
470 else
471 {
474 }
475 }
476 else
477 {
478 if (IsForwards)
479 {
482 }
483 else
484 {
487 }
488 }
489
490 int nrowsY = blkmatY->GetRows();
491 int ncolsY = blkmatY->GetColumns();
492
493 Array<OneD, NekDouble> sortedinarrayY(ncolsY);
494 Array<OneD, NekDouble> sortedoutarrayY(nrowsY);
495
496 int nrowsZ = blkmatZ->GetRows();
497 int ncolsZ = blkmatZ->GetColumns();
498
499 Array<OneD, NekDouble> sortedinarrayZ(ncolsZ);
500 Array<OneD, NekDouble> sortedoutarrayZ(nrowsZ);
501
502 NekVector<NekDouble> inY(ncolsY, sortedinarrayY, eWrapper);
503 NekVector<NekDouble> outY(nrowsY, sortedoutarrayY, eWrapper);
504
505 NekVector<NekDouble> inZ(ncolsZ, sortedinarrayZ, eWrapper);
506 NekVector<NekDouble> outZ(nrowsZ, sortedoutarrayZ, eWrapper);
507
508 m_transposition->Transpose(num_dofs, inarray, sortedinarrayY,
509 !IsForwards, LibUtilities::eXtoYZ);
510
511 outY = (*blkmatY) * inY;
512
513 m_transposition->Transpose(sortedoutarrayY.size(), sortedoutarrayY,
514 sortedinarrayZ, false,
516
517 outZ = (*blkmatZ) * inZ;
518
519 m_transposition->Transpose(sortedoutarrayZ.size(), sortedoutarrayZ,
520 sortedoutarrayY, false,
522
523 m_transposition->Transpose(sortedoutarrayY.size(), sortedoutarrayY,
524 outarray, false, LibUtilities::eYZtoX);
525 }
526}
DNekBlkMatSharedPtr GetHomogeneous2DBlockMatrix(Homogeneous2DMatType mattype) const

References Nektar::MultiRegions::eBackwardsCoeffSpaceY1D, Nektar::MultiRegions::eBackwardsCoeffSpaceZ1D, Nektar::MultiRegions::eBackwardsPhysSpaceY1D, Nektar::MultiRegions::eBackwardsPhysSpaceZ1D, Nektar::MultiRegions::eForwardsCoeffSpaceY1D, Nektar::MultiRegions::eForwardsCoeffSpaceZ1D, Nektar::MultiRegions::eForwardsPhysSpaceY1D, Nektar::MultiRegions::eForwardsPhysSpaceZ1D, Nektar::eWrapper, Nektar::LibUtilities::eXtoYZ, Nektar::LibUtilities::eYZtoX, Nektar::LibUtilities::eYZtoZY, Nektar::LibUtilities::eZYtoYZ, GetHomogeneous2DBlockMatrix(), m_FFT_y, m_FFT_z, m_lines, Nektar::MultiRegions::ExpList::m_npoints, m_ny, m_nz, m_tmpIN, m_tmpOUT, m_transposition, m_useFFT, and CellMLToNektar.cellml_metadata::p.

Referenced by v_HomogeneousBwdTrans(), and v_HomogeneousFwdTrans().

◆ PhysDeriv() [1/2]

void Nektar::MultiRegions::ExpListHomogeneous2D::PhysDeriv ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d0,
Array< OneD, NekDouble > &  out_d1,
Array< OneD, NekDouble > &  out_d2 
)

Definition at line 1057 of file ExpListHomogeneous2D.cpp.

1061{
1062 v_PhysDeriv(inarray, out_d0, out_d1, out_d2);
1063}
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

References v_PhysDeriv().

◆ PhysDeriv() [2/2]

void Nektar::MultiRegions::ExpListHomogeneous2D::PhysDeriv ( Direction  edir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d 
)

Definition at line 1065 of file ExpListHomogeneous2D.cpp.

1068{
1069 // convert int into enum
1070 v_PhysDeriv(edir, inarray, out_d);
1071}

References v_PhysDeriv().

◆ SetPaddingBase()

void Nektar::MultiRegions::ExpListHomogeneous2D::SetPaddingBase ( void  )

Definition at line 1073 of file ExpListHomogeneous2D.cpp.

1074{
1075 NekDouble size_y = 1.5 * m_ny;
1076 NekDouble size_z = 1.5 * m_nz;
1077 m_padsize_y = int(size_y);
1078 m_padsize_z = int(size_z);
1079
1080 const LibUtilities::PointsKey Ppad_y(m_padsize_y,
1082 const LibUtilities::BasisKey Bpad_y(LibUtilities::eFourier, m_padsize_y,
1083 Ppad_y);
1084
1085 const LibUtilities::PointsKey Ppad_z(m_padsize_z,
1087 const LibUtilities::BasisKey Bpad_z(LibUtilities::eFourier, m_padsize_z,
1088 Ppad_z);
1089
1092
1093 StdRegions::StdQuadExp StdQuad(m_paddingBasis_y->GetBasisKey(),
1094 m_paddingBasis_z->GetBasisKey());
1095
1096 StdRegions::StdMatrixKey matkey1(StdRegions::eFwdTrans,
1097 StdQuad.DetShapeType(), StdQuad);
1098 StdRegions::StdMatrixKey matkey2(StdRegions::eBwdTrans,
1099 StdQuad.DetShapeType(), StdQuad);
1100
1101 MatFwdPAD = StdQuad.GetStdMatrix(matkey1);
1102 MatBwdPAD = StdQuad.GetStdMatrix(matkey2);
1103}
LibUtilities::BasisSharedPtr m_paddingBasis_z
Base expansion in z direction.
LibUtilities::BasisSharedPtr m_paddingBasis_y
Base expansion in y direction.
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:74
@ eFourier
Fourier Expansion .
Definition: BasisType.h:55
double NekDouble

References Nektar::LibUtilities::BasisManager(), Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::StdRegions::eBwdTrans, Nektar::LibUtilities::eFourier, Nektar::LibUtilities::eFourierEvenlySpaced, Nektar::StdRegions::eFwdTrans, Nektar::StdRegions::StdExpansion::GetStdMatrix(), m_ny, m_nz, m_paddingBasis_y, m_paddingBasis_z, m_padsize_y, m_padsize_z, MatBwdPAD, and MatFwdPAD.

Referenced by ExpListHomogeneous2D().

◆ v_AppendFieldData() [1/2]

void Nektar::MultiRegions::ExpListHomogeneous2D::v_AppendFieldData ( LibUtilities::FieldDefinitionsSharedPtr fielddef,
std::vector< NekDouble > &  fielddata 
)
overrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 741 of file ExpListHomogeneous2D.cpp.

744{
745 v_AppendFieldData(fielddef, fielddata, m_coeffs);
746}
void v_AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata) override
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:1080

References Nektar::MultiRegions::ExpList::m_coeffs, and v_AppendFieldData().

Referenced by v_AppendFieldData().

◆ v_AppendFieldData() [2/2]

void Nektar::MultiRegions::ExpListHomogeneous2D::v_AppendFieldData ( LibUtilities::FieldDefinitionsSharedPtr fielddef,
std::vector< NekDouble > &  fielddata,
Array< OneD, NekDouble > &  coeffs 
)
overrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 707 of file ExpListHomogeneous2D.cpp.

710{
711 int i, k;
712
713 int NumMod_y = m_homogeneousBasis_y->GetNumModes();
714 int NumMod_z = m_homogeneousBasis_z->GetNumModes();
715
716 int ncoeffs_per_line = m_lines[0]->GetNcoeffs();
717
718 // Determine mapping from element ids to location in
719 // expansion list
720 map<int, int> ElmtID_to_ExpID;
721 for (i = 0; i < m_lines[0]->GetExpSize(); ++i)
722 {
723 ElmtID_to_ExpID[(*m_exp)[i]->GetGeom()->GetGlobalID()] = i;
724 }
725
726 for (i = 0; i < fielddef->m_elementIDs.size(); ++i)
727 {
728 int eid = ElmtID_to_ExpID[fielddef->m_elementIDs[i]];
729 int datalen = (*m_exp)[eid]->GetNcoeffs();
730
731 for (k = 0; k < (NumMod_y * NumMod_z); ++k)
732 {
733 fielddata.insert(
734 fielddata.end(),
735 &coeffs[m_coeff_offset[eid] + k * ncoeffs_per_line],
736 &coeffs[m_coeff_offset[eid] + k * ncoeffs_per_line] + datalen);
737 }
738 }
739}
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:1124

References Nektar::MultiRegions::ExpList::m_coeff_offset, m_homogeneousBasis_y, m_homogeneousBasis_z, and m_lines.

◆ v_BwdTrans()

void Nektar::MultiRegions::ExpListHomogeneous2D::v_BwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Given the elemental coefficients \(\hat{u}_n^e\) of an expansion, this function evaluates the spectral/hp expansion \(u^{\delta}(\boldsymbol{x})\) at the quadrature points \(\boldsymbol{x}_i\). The operation is evaluated locally by the elemental function StdRegions::StdExpansion::BwdTrans.

Parameters
inarrayAn array of size \(N_{\mathrm{eof}}\) containing the local coefficients \(\hat{u}_n^e\).
outarrayThe resulting physical values at the quadrature points \(u^{\delta}(\boldsymbol{x}_i)\) will be stored in this array of size \(Q_{\mathrm{tot}}\).

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 357 of file ExpListHomogeneous2D.cpp.

360{
361 size_t cnt = 0, cnt1 = 0;
362 Array<OneD, NekDouble> tmparray;
363 size_t nlines = m_lines.size();
364
365 for (size_t n = 0; n < nlines; ++n)
366 {
367 m_lines[n]->BwdTrans(inarray + cnt, tmparray = outarray + cnt1);
368 cnt += m_lines[n]->GetNcoeffs();
369 cnt1 += m_lines[n]->GetTotPoints();
370 }
371 if (!m_WaveSpace)
372 {
373 HomogeneousBwdTrans(cnt1, outarray, outarray);
374 }
375}
void HomogeneousBwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:1857

References Nektar::MultiRegions::ExpList::HomogeneousBwdTrans(), m_lines, and Nektar::MultiRegions::ExpList::m_WaveSpace.

◆ v_DealiasedDotProd()

void Nektar::MultiRegions::ExpListHomogeneous2D::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 
)
overrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 269 of file ExpListHomogeneous2D.cpp.

273{
274 // TODO Proper implementation of this
275 size_t ndim = inarray1.size();
276 ASSERTL1(inarray2.size() % ndim == 0,
277 "Wrong dimensions for DealiasedDotProd.");
278
279 size_t nvec = inarray2.size() % ndim;
280
281 Array<OneD, NekDouble> out(npts);
282 for (size_t i = 0; i < nvec; i++)
283 {
284 Vmath::Zero(npts, outarray[i], 1);
285 for (size_t j = 0; j < ndim; j++)
286 {
287 DealiasedProd(npts, inarray1[j], inarray2[i * ndim + j], out);
288 Vmath::Vadd(npts, outarray[i], 1, out, 1, outarray[i], 1);
289 }
290 }
291}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
void DealiasedProd(const int num_dofs, const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1866
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273

References ASSERTL1, Nektar::MultiRegions::ExpList::DealiasedProd(), Vmath::Vadd(), and Vmath::Zero().

◆ v_DealiasedProd()

void Nektar::MultiRegions::ExpListHomogeneous2D::v_DealiasedProd ( const int  num_dofs,
const Array< OneD, NekDouble > &  inarray1,
const Array< OneD, NekDouble > &  inarray2,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 162 of file ExpListHomogeneous2D.cpp.

165{
166 // int npoints = outarray.size(); // number of total physical points
167 int nlines =
168 m_lines.size(); // number of lines == number of Fourier modes = number
169 // of Fourier coeff = number of points per slab
170 int nslabs = npoints /
171 nlines; // number of slabs = numebr of physical points per line
172
173 Array<OneD, NekDouble> V1(npoints);
174 Array<OneD, NekDouble> V2(npoints);
175 Array<OneD, NekDouble> V1V2(npoints);
176 Array<OneD, NekDouble> ShufV1(npoints);
177 Array<OneD, NekDouble> ShufV2(npoints);
178 Array<OneD, NekDouble> ShufV1V2(npoints);
179
180 if (m_WaveSpace)
181 {
182 V1 = inarray1;
183 V2 = inarray2;
184 }
185 else
186 {
187 HomogeneousFwdTrans(npoints, inarray1, V1);
188 HomogeneousFwdTrans(npoints, inarray2, V2);
189 }
190
191 m_transposition->Transpose(npoints, V1, ShufV1, false,
193 m_transposition->Transpose(npoints, V2, ShufV2, false,
195
196 Array<OneD, NekDouble> PadV1_slab_coeff(m_padsize_y * m_padsize_z, 0.0);
197 Array<OneD, NekDouble> PadV2_slab_coeff(m_padsize_y * m_padsize_z, 0.0);
198 Array<OneD, NekDouble> PadRe_slab_coeff(m_padsize_y * m_padsize_z, 0.0);
199
200 Array<OneD, NekDouble> PadV1_slab_phys(m_padsize_y * m_padsize_z, 0.0);
201 Array<OneD, NekDouble> PadV2_slab_phys(m_padsize_y * m_padsize_z, 0.0);
202 Array<OneD, NekDouble> PadRe_slab_phys(m_padsize_y * m_padsize_z, 0.0);
203
204 NekVector<NekDouble> PadIN_V1(m_padsize_y * m_padsize_z, PadV1_slab_coeff,
205 eWrapper);
206 NekVector<NekDouble> PadOUT_V1(m_padsize_y * m_padsize_z, PadV1_slab_phys,
207 eWrapper);
208
209 NekVector<NekDouble> PadIN_V2(m_padsize_y * m_padsize_z, PadV2_slab_coeff,
210 eWrapper);
211 NekVector<NekDouble> PadOUT_V2(m_padsize_y * m_padsize_z, PadV2_slab_phys,
212 eWrapper);
213
214 NekVector<NekDouble> PadIN_Re(m_padsize_y * m_padsize_z, PadRe_slab_phys,
215 eWrapper);
216 NekVector<NekDouble> PadOUT_Re(m_padsize_y * m_padsize_z, PadRe_slab_coeff,
217 eWrapper);
218
219 // Looping on the slabs
220 for (int j = 0; j < nslabs; j++)
221 {
222 // Copying the j-th slab of size N*M into a bigger slab of lenght 2*N*M
223 // We are in Fourier space
224 for (int i = 0; i < m_nz; i++)
225 {
226 Vmath::Vcopy(m_ny, &(ShufV1[i * m_ny + j * nlines]), 1,
227 &(PadV1_slab_coeff[i * 2 * m_ny]), 1);
228 Vmath::Vcopy(m_ny, &(ShufV2[i * m_ny + j * nlines]), 1,
229 &(PadV2_slab_coeff[i * 2 * m_ny]), 1);
230 }
231
232 // Moving to physical space using the padded system
233 PadOUT_V1 = (*MatBwdPAD) * PadIN_V1;
234 PadOUT_V2 = (*MatBwdPAD) * PadIN_V2;
235
236 // Perfroming the vectors multiplication in physical
237 // space on the padded system
238 Vmath::Vmul(m_padsize_y * m_padsize_z, PadV1_slab_phys, 1,
239 PadV2_slab_phys, 1, PadRe_slab_phys, 1);
240
241 // Moving back the result (V1*V2)_phys in Fourier
242 // space, padded system
243 PadOUT_Re = (*MatFwdPAD) * PadIN_Re;
244
245 // Copying the first half of the padded pencil in the
246 // full vector (Fourier space)
247 for (int i = 0; i < m_nz; i++)
248 {
249 Vmath::Vcopy(m_ny, &(PadRe_slab_coeff[i * 2 * m_ny]), 1,
250 &(ShufV1V2[i * m_ny + j * nlines]), 1);
251 }
252 }
253
254 if (m_WaveSpace)
255 {
256 m_transposition->Transpose(npoints, ShufV1V2, outarray, false,
258 }
259 else
260 {
261 m_transposition->Transpose(npoints, ShufV1V2, V1V2, false,
263
264 // Moving the results in physical space for the output
265 HomogeneousBwdTrans(npoints, V1V2, outarray);
266 }
267}
void HomogeneousFwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:1848
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.hpp:72
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References Nektar::eWrapper, Nektar::LibUtilities::eXtoYZ, Nektar::LibUtilities::eYZtoX, Nektar::MultiRegions::ExpList::HomogeneousBwdTrans(), Nektar::MultiRegions::ExpList::HomogeneousFwdTrans(), m_lines, m_ny, m_nz, m_padsize_y, m_padsize_z, m_transposition, Nektar::MultiRegions::ExpList::m_WaveSpace, Vmath::Vcopy(), and Vmath::Vmul().

◆ v_ExtractDataToCoeffs()

void Nektar::MultiRegions::ExpListHomogeneous2D::v_ExtractDataToCoeffs ( LibUtilities::FieldDefinitionsSharedPtr fielddef,
std::vector< NekDouble > &  fielddata,
std::string &  field,
Array< OneD, NekDouble > &  coeffs,
std::unordered_map< int, int >  zIdToPlane 
)
overrideprotectedvirtual

Extract data from raw field data into expansion list.

Parameters
fielddefField definitions.
fielddataData for associated field.
fieldField variable name.
coeffsResulting coefficient array.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 748 of file ExpListHomogeneous2D.cpp.

753{
754 int i, k;
755 int offset = 0;
756 int datalen = fielddata.size() / fielddef->m_fields.size();
757 int ncoeffs_per_line = m_lines[0]->GetNcoeffs();
758 int NumMod_y = m_homogeneousBasis_y->GetNumModes();
759 int NumMod_z = m_homogeneousBasis_z->GetNumModes();
760 // Find data location according to field definition
761 for (i = 0; i < fielddef->m_fields.size(); ++i)
762 {
763 if (fielddef->m_fields[i] == field)
764 {
765 break;
766 }
767 offset += datalen;
768 }
769 ASSERTL0(i != fielddef->m_fields.size(), "Field not found in data file");
770 // Determine mapping from element ids to location in
771 // expansion list
772 map<int, int> ElmtID_to_ExpID;
773 for (i = 0; i < m_lines[0]->GetExpSize(); ++i)
774 {
775 ElmtID_to_ExpID[(*m_exp)[i]->GetGeom()->GetGlobalID()] = i;
776 }
777 for (i = 0; i < fielddef->m_elementIDs.size(); ++i)
778 {
779 int eid = ElmtID_to_ExpID[fielddef->m_elementIDs[i]];
780 int datalen = (*m_exp)[eid]->GetNcoeffs();
781 for (k = 0; k < (NumMod_y * NumMod_z); ++k)
782 {
783 Vmath::Vcopy(datalen, &fielddata[offset], 1,
784 &coeffs[m_coeff_offset[eid] + k * ncoeffs_per_line],
785 1);
786 offset += datalen;
787 }
788 }
789}

References ASSERTL0, Nektar::MultiRegions::ExpList::m_coeff_offset, m_homogeneousBasis_y, m_homogeneousBasis_z, m_lines, and Vmath::Vcopy().

◆ v_FwdTrans()

void Nektar::MultiRegions::ExpListHomogeneous2D::v_FwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 293 of file ExpListHomogeneous2D.cpp.

296{
297 size_t cnt = 0, cnt1 = 0;
298 Array<OneD, NekDouble> tmparray;
299 size_t nlines = m_lines.size();
300
301 for (size_t n = 0; n < nlines; ++n)
302 {
303 m_lines[n]->FwdTrans(inarray + cnt, tmparray = outarray + cnt1);
304 cnt += m_lines[n]->GetTotPoints();
305 cnt1 += m_lines[n]->GetNcoeffs();
306 }
307 if (!m_WaveSpace)
308 {
309 HomogeneousFwdTrans(cnt1, outarray, outarray);
310 }
311}

References Nektar::MultiRegions::ExpList::HomogeneousFwdTrans(), m_lines, and Nektar::MultiRegions::ExpList::m_WaveSpace.

◆ v_FwdTransBndConstrained()

void Nektar::MultiRegions::ExpListHomogeneous2D::v_FwdTransBndConstrained ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 335 of file ExpListHomogeneous2D.cpp.

338{
339 int cnt = 0, cnt1 = 0;
340 Array<OneD, NekDouble> tmparray;
341 int nlines = m_lines.size();
342
343 for (int n = 0; n < nlines; ++n)
344 {
345 m_lines[n]->FwdTransBndConstrained(inarray + cnt,
346 tmparray = outarray + cnt1);
347
348 cnt += m_lines[n]->GetTotPoints();
349 cnt1 += m_lines[n]->GetNcoeffs();
350 }
351 if (!m_WaveSpace)
352 {
353 HomogeneousFwdTrans(cnt1, outarray, outarray);
354 }
355}

References Nektar::MultiRegions::ExpList::HomogeneousFwdTrans(), m_lines, and Nektar::MultiRegions::ExpList::m_WaveSpace.

◆ v_FwdTransLocalElmt()

void Nektar::MultiRegions::ExpListHomogeneous2D::v_FwdTransLocalElmt ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

Given a function \(u(\boldsymbol{x})\) defined at the quadrature points, this function determines the transformed elemental coefficients \(\hat{u}_n^e\) employing a discrete elemental Galerkin projection from physical space to coefficient space. For each element, the operation is evaluated locally by the function StdRegions::StdExpansion::IproductWRTBase followed by a call to #MultiRegions#MultiplyByElmtInvMass.

Parameters
inarrayAn array of size \(Q_{\mathrm{tot}}\) containing the values of the function \(f(\boldsymbol{x})\) at the quadrature points \(\boldsymbol{x}_i\).
outarrayThe resulting coefficients \(\hat{u}_n^e\) will be stored in this array of size \(N_{\mathrm{eof}}\).

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 313 of file ExpListHomogeneous2D.cpp.

316{
317 size_t cnt = 0, cnt1 = 0;
318 Array<OneD, NekDouble> tmparray;
319 size_t nlines = m_lines.size();
320
321 for (size_t n = 0; n < nlines; ++n)
322 {
323 m_lines[n]->FwdTransLocalElmt(inarray + cnt,
324 tmparray = outarray + cnt1);
325
326 cnt += m_lines[n]->GetTotPoints();
327 cnt1 += m_lines[n]->GetNcoeffs();
328 }
329 if (!m_WaveSpace)
330 {
331 HomogeneousFwdTrans(cnt1, outarray, outarray);
332 }
333}

References Nektar::MultiRegions::ExpList::HomogeneousFwdTrans(), m_lines, and Nektar::MultiRegions::ExpList::m_WaveSpace.

◆ v_GetFieldDefinitions() [1/2]

void Nektar::MultiRegions::ExpListHomogeneous2D::v_GetFieldDefinitions ( std::vector< LibUtilities::FieldDefinitionsSharedPtr > &  fielddef)
overrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 674 of file ExpListHomogeneous2D.cpp.

676{
677 // Set up Homogeneous length details.
678 Array<OneD, LibUtilities::BasisSharedPtr> HomoBasis(2);
679 HomoBasis[0] = m_homogeneousBasis_y;
680 HomoBasis[1] = m_homogeneousBasis_z;
681 std::vector<NekDouble> HomoLen(2);
682 HomoLen[0] = m_lhom_y;
683 HomoLen[1] = m_lhom_z;
684
685 int nhom_modes_y = m_homogeneousBasis_y->GetNumModes();
686 int nhom_modes_z = m_homogeneousBasis_z->GetNumModes();
687
688 std::vector<unsigned int> sIDs = LibUtilities::NullUnsignedIntVector;
689
690 std::vector<unsigned int> yIDs;
691 std::vector<unsigned int> zIDs;
692
693 for (int n = 0; n < nhom_modes_z; ++n)
694 {
695 for (int m = 0; m < nhom_modes_y; ++m)
696 {
697 zIDs.push_back(n);
698 yIDs.push_back(m);
699 }
700 }
701
702 // enforce NumHomoDir == 1 by direct call
703 m_lines[0]->GeneralGetFieldDefinitions(fielddef, 2, HomoBasis, HomoLen,
704 false, sIDs, zIDs, yIDs);
705}
static std::vector< unsigned int > NullUnsignedIntVector

References m_homogeneousBasis_y, m_homogeneousBasis_z, m_lhom_y, m_lhom_z, m_lines, and Nektar::LibUtilities::NullUnsignedIntVector.

◆ v_GetFieldDefinitions() [2/2]

std::vector< LibUtilities::FieldDefinitionsSharedPtr > Nektar::MultiRegions::ExpListHomogeneous2D::v_GetFieldDefinitions ( void  )
overrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 639 of file ExpListHomogeneous2D.cpp.

641{
642 std::vector<LibUtilities::FieldDefinitionsSharedPtr> returnval;
643 // Set up Homogeneous length details.
644 Array<OneD, LibUtilities::BasisSharedPtr> HomoBasis(2);
645 HomoBasis[0] = m_homogeneousBasis_y;
646 HomoBasis[1] = m_homogeneousBasis_z;
647
648 std::vector<NekDouble> HomoLen(2);
649 HomoLen[0] = m_lhom_y;
650 HomoLen[1] = m_lhom_z;
651
652 int nhom_modes_y = m_homogeneousBasis_y->GetNumModes();
653 int nhom_modes_z = m_homogeneousBasis_z->GetNumModes();
654
655 std::vector<unsigned int> sIDs = LibUtilities::NullUnsignedIntVector;
656
657 std::vector<unsigned int> yIDs;
658 std::vector<unsigned int> zIDs;
659
660 for (int n = 0; n < nhom_modes_z; ++n)
661 {
662 for (int m = 0; m < nhom_modes_y; ++m)
663 {
664 zIDs.push_back(n);
665 yIDs.push_back(m);
666 }
667 }
668
669 m_lines[0]->GeneralGetFieldDefinitions(returnval, 2, HomoBasis, HomoLen,
670 false, sIDs, zIDs, yIDs);
671 return returnval;
672}

References m_homogeneousBasis_y, m_homogeneousBasis_z, m_lhom_y, m_lhom_z, m_lines, and Nektar::LibUtilities::NullUnsignedIntVector.

◆ v_GetNumElmts()

size_t Nektar::MultiRegions::ExpListHomogeneous2D::v_GetNumElmts ( void  )
inlineoverrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 152 of file ExpListHomogeneous2D.h.

153 {
154 return m_lines[0]->GetExpSize();
155 }

References m_lines.

◆ v_HomogeneousBwdTrans()

void Nektar::MultiRegions::ExpListHomogeneous2D::v_HomogeneousBwdTrans ( const int  npts,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
bool  Shuff = true,
bool  UnShuff = true 
)
overrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 154 of file ExpListHomogeneous2D.cpp.

157{
158 // Backwards trans
159 Homogeneous2DTrans(npts, inarray, outarray, false, Shuff, UnShuff);
160}
void Homogeneous2DTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool IsForwards, bool Shuff=true, bool UnShuff=true)

References Homogeneous2DTrans().

◆ v_HomogeneousFwdTrans()

void Nektar::MultiRegions::ExpListHomogeneous2D::v_HomogeneousFwdTrans ( const int  npts,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
bool  Shuff = true,
bool  UnShuff = true 
)
overrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 146 of file ExpListHomogeneous2D.cpp.

149{
150 // Forwards trans
151 Homogeneous2DTrans(npts, inarray, outarray, true, Shuff, UnShuff);
152}

References Homogeneous2DTrans().

◆ v_IProductWRTBase()

void Nektar::MultiRegions::ExpListHomogeneous2D::v_IProductWRTBase ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
overrideprotectedvirtual

The operation is evaluated locally for every element by the function StdRegions::StdExpansion::IProductWRTBase.

Parameters
inarrayAn array of size \(Q_{\mathrm{tot}}\) containing the values of the function \(f(\boldsymbol{x})\) at the quadrature points \(\boldsymbol{x}_i\).
outarrayAn array of size \(N_{\mathrm{eof}}\) used to store the result.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 377 of file ExpListHomogeneous2D.cpp.

380{
381 size_t cnt = 0, cnt1 = 0;
382 Array<OneD, NekDouble> tmparray;
383 size_t nlines = m_lines.size();
384
385 for (size_t n = 0; n < nlines; ++n)
386 {
387 m_lines[n]->IProductWRTBase(inarray + cnt, tmparray = outarray + cnt1);
388
389 cnt += m_lines[n]->GetNcoeffs();
390 cnt1 += m_lines[n]->GetTotPoints();
391 }
392}

References m_lines.

◆ v_PhysDeriv() [1/2]

void Nektar::MultiRegions::ExpListHomogeneous2D::v_PhysDeriv ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d0,
Array< OneD, NekDouble > &  out_d1,
Array< OneD, NekDouble > &  out_d2 
)
overrideprotectedvirtual

Given a function \(f(\boldsymbol{x})\) evaluated at the quadrature points, this function calculates the derivatives \(\frac{d}{dx_1}\), \(\frac{d}{dx_2}\) and \(\frac{d}{dx_3}\) of the function \(f(\boldsymbol{x})\) at the same quadrature points. The local distribution of the quadrature points allows an elemental evaluation of the derivative. This is done by a call to the function StdRegions::StdExpansion::PhysDeriv.

Parameters
inarrayAn array of size \(Q_{\mathrm{tot}}\) containing the values of the function \(f(\boldsymbol{x})\) at the quadrature points \(\boldsymbol{x}_i\).
out_d0The discrete evaluation of the derivative \(\frac{d}{dx_1}\) will be stored in this array of size \(Q_{\mathrm{tot}}\).
out_d1The discrete evaluation of the derivative \(\frac{d}{dx_2}\) will be stored in this array of size \(Q_{\mathrm{tot}}\). Note that if no memory is allocated for out_d1, the derivative \(\frac{d}{dx_2}\) will not be calculated.
out_d2The discrete evaluation of the derivative \(\frac{d}{dx_3}\) will be stored in this array of size \(Q_{\mathrm{tot}}\). Note that if no memory is allocated for out_d2, the derivative \(\frac{d}{dx_3}\) will not be calculated.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 816 of file ExpListHomogeneous2D.cpp.

820{
821 int nyzlines = m_lines.size(); // number of Fourier points in the Fourier
822 // directions (nF_pts)
823 int n_points_line = m_npoints / nyzlines; // number of points per line
824
825 Array<OneD, NekDouble> temparray(m_npoints);
826 Array<OneD, NekDouble> temparray1(m_npoints);
827 Array<OneD, NekDouble> temparray2(m_npoints);
828 Array<OneD, NekDouble> tmp1;
829 Array<OneD, NekDouble> tmp2;
830 Array<OneD, NekDouble> tmp3;
831
832 for (int i = 0; i < nyzlines; i++)
833 {
834 m_lines[i]->PhysDeriv(tmp1 = inarray + i * n_points_line,
835 tmp2 = out_d0 + i * n_points_line);
836 }
837
838 if (m_homogeneousBasis_y->GetBasisType() == LibUtilities::eFourier &&
840 {
841 if (m_WaveSpace)
842 {
843 temparray = inarray;
844 }
845 else
846 {
847 HomogeneousFwdTrans(m_npoints, inarray, temparray);
848 }
849 NekDouble sign = -1.0;
851
852 // along y
853 for (int i = 0; i < m_ny; i++)
854 {
855 beta = -sign * 2 * M_PI * (i / 2) / m_lhom_y;
856
857 for (int j = 0; j < m_nz; j++)
858 {
859 Vmath::Smul(n_points_line, beta,
860 tmp1 = temparray + n_points_line * (i + j * m_ny),
861 1,
862 tmp2 = temparray1 +
863 n_points_line * ((i - int(sign)) + j * m_ny),
864 1);
865 }
866
867 sign = -1.0 * sign;
868 }
869
870 // along z
871 sign = -1.0;
872 for (int i = 0; i < m_nz; i++)
873 {
874 beta = -sign * 2 * M_PI * (i / 2) / m_lhom_z;
876 m_ny * n_points_line, beta,
877 tmp1 = temparray + i * m_ny * n_points_line, 1,
878 tmp2 = temparray2 + (i - int(sign)) * m_ny * n_points_line, 1);
879 sign = -1.0 * sign;
880 }
881 if (m_WaveSpace)
882 {
883 out_d1 = temparray1;
884 out_d2 = temparray2;
885 }
886 else
887 {
888 HomogeneousBwdTrans(m_npoints, temparray1, out_d1);
889 HomogeneousBwdTrans(m_npoints, temparray2, out_d2);
890 }
891 }
892 else
893 {
894 if (m_WaveSpace)
895 {
896 ASSERTL0(false, "Semi-phyisical time-stepping not implemented yet "
897 "for non-Fourier basis")
898 }
899 else
900 {
901 StdRegions::StdQuadExp StdQuad(m_homogeneousBasis_y->GetBasisKey(),
902 m_homogeneousBasis_z->GetBasisKey());
903
904 m_transposition->Transpose(m_npoints, inarray, temparray, false,
906
907 for (int i = 0; i < n_points_line; i++)
908 {
909 StdQuad.PhysDeriv(tmp1 = temparray + i * nyzlines,
910 tmp2 = temparray1 + i * nyzlines,
911 tmp3 = temparray2 + i * nyzlines);
912 }
913
914 m_transposition->Transpose(m_npoints, temparray1, out_d1, false,
916 m_transposition->Transpose(m_npoints, temparray2, out_d2, false,
918 Vmath::Smul(m_npoints, 2.0 / m_lhom_y, out_d1, 1, out_d1, 1);
919 Vmath::Smul(m_npoints, 2.0 / m_lhom_z, out_d2, 1, out_d2, 1);
920 }
921 }
922}
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:47
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:59
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

References ASSERTL0, Nektar::LibUtilities::beta, Nektar::LibUtilities::eFourier, Nektar::LibUtilities::eXtoYZ, Nektar::LibUtilities::eYZtoX, Nektar::MultiRegions::ExpList::HomogeneousBwdTrans(), Nektar::MultiRegions::ExpList::HomogeneousFwdTrans(), m_homogeneousBasis_y, m_homogeneousBasis_z, m_lhom_y, m_lhom_z, m_lines, Nektar::MultiRegions::ExpList::m_npoints, m_ny, m_nz, m_transposition, Nektar::MultiRegions::ExpList::m_WaveSpace, Nektar::StdRegions::StdExpansion::PhysDeriv(), sign, and Vmath::Smul().

Referenced by PhysDeriv().

◆ v_PhysDeriv() [2/2]

void Nektar::MultiRegions::ExpListHomogeneous2D::v_PhysDeriv ( Direction  edir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d 
)
overrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 924 of file ExpListHomogeneous2D.cpp.

928{
929 int nyzlines = m_lines.size(); // number of Fourier points in the Fourier
930 // directions (nF_pts)
931 int n_points_line = m_npoints / nyzlines; // number of points per line
932 // convert enum into int
933 int dir = (int)edir;
934
935 Array<OneD, NekDouble> temparray(m_npoints);
936 Array<OneD, NekDouble> temparray1(m_npoints);
937 Array<OneD, NekDouble> temparray2(m_npoints);
938 Array<OneD, NekDouble> tmp1;
939 Array<OneD, NekDouble> tmp2;
940 Array<OneD, NekDouble> tmp3;
941
942 if (dir < 1)
943 {
944 for (int i = 0; i < nyzlines; i++)
945 {
946 m_lines[i]->PhysDeriv(tmp1 = inarray + i * n_points_line,
947 tmp2 = out_d + i * n_points_line);
948 }
949 }
950 else
951 {
952 if (m_homogeneousBasis_y->GetBasisType() == LibUtilities::eFourier &&
954 {
955 if (m_WaveSpace)
956 {
957 temparray = inarray;
958 }
959 else
960 {
961 HomogeneousFwdTrans(m_npoints, inarray, temparray);
962 }
963 NekDouble sign = -1.0;
965
966 if (dir == 1)
967 {
968 // along y
969 for (int i = 0; i < m_ny; i++)
970 {
971 beta = -sign * 2 * M_PI * (i / 2) / m_lhom_y;
972
973 for (int j = 0; j < m_nz; j++)
974 {
976 n_points_line, beta,
977 tmp1 = temparray + n_points_line * (i + j * m_ny),
978 1,
979 tmp2 = temparray1 +
980 n_points_line * ((i - int(sign)) + j * m_ny),
981 1);
982 }
983 sign = -1.0 * sign;
984 }
985 if (m_WaveSpace)
986 {
987 out_d = temparray1;
988 }
989 else
990 {
991 HomogeneousBwdTrans(m_npoints, temparray1, out_d);
992 }
993 }
994 else
995 {
996 // along z
997 for (int i = 0; i < m_nz; i++)
998 {
999 beta = -sign * 2 * M_PI * (i / 2) / m_lhom_z;
1000 Vmath::Smul(m_ny * n_points_line, beta,
1001 tmp1 = temparray + i * m_ny * n_points_line, 1,
1002 tmp2 = temparray2 +
1003 (i - int(sign)) * m_ny * n_points_line,
1004 1);
1005 sign = -1.0 * sign;
1006 }
1007 if (m_WaveSpace)
1008 {
1009 out_d = temparray2;
1010 }
1011 else
1012 {
1013 HomogeneousBwdTrans(m_npoints, temparray2, out_d);
1014 }
1015 }
1016 }
1017 else
1018 {
1019 if (m_WaveSpace)
1020 {
1021 ASSERTL0(false, "Semi-phyisical time-stepping not implemented "
1022 "yet for non-Fourier basis")
1023 }
1024 else
1025 {
1026 StdRegions::StdQuadExp StdQuad(
1027 m_homogeneousBasis_y->GetBasisKey(),
1028 m_homogeneousBasis_z->GetBasisKey());
1029
1030 m_transposition->Transpose(m_npoints, inarray, temparray, false,
1032
1033 for (int i = 0; i < n_points_line; i++)
1034 {
1035 StdQuad.PhysDeriv(tmp1 = temparray + i * nyzlines,
1036 tmp2 = temparray1 + i * nyzlines,
1037 tmp3 = temparray2 + i * nyzlines);
1038 }
1039
1040 if (dir == 1)
1041 {
1042 m_transposition->Transpose(m_npoints, temparray1, out_d,
1043 false, LibUtilities::eYZtoX);
1044 Vmath::Smul(m_npoints, 2.0 / m_lhom_y, out_d, 1, out_d, 1);
1045 }
1046 else
1047 {
1048 m_transposition->Transpose(m_npoints, temparray2, out_d,
1049 false, LibUtilities::eYZtoX);
1050 Vmath::Smul(m_npoints, 2.0 / m_lhom_z, out_d, 1, out_d, 1);
1051 }
1052 }
1053 }
1054 }
1055}

References ASSERTL0, Nektar::LibUtilities::beta, Nektar::LibUtilities::eFourier, Nektar::LibUtilities::eXtoYZ, Nektar::LibUtilities::eYZtoX, Nektar::MultiRegions::ExpList::HomogeneousBwdTrans(), Nektar::MultiRegions::ExpList::HomogeneousFwdTrans(), m_homogeneousBasis_y, m_homogeneousBasis_z, m_lhom_y, m_lhom_z, m_lines, Nektar::MultiRegions::ExpList::m_npoints, m_ny, m_nz, m_transposition, Nektar::MultiRegions::ExpList::m_WaveSpace, Nektar::StdRegions::StdExpansion::PhysDeriv(), sign, and Vmath::Smul().

◆ v_WriteVtkPieceData()

void Nektar::MultiRegions::ExpListHomogeneous2D::v_WriteVtkPieceData ( std::ostream &  outfile,
int  expansion,
std::string  var 
)
overrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 790 of file ExpListHomogeneous2D.cpp.

792{
793 int i;
794 int nq = (*m_exp)[expansion]->GetTotPoints();
795 int npoints_per_line = m_lines[0]->GetTotPoints();
796
797 // printing the fields of that zone
798 outfile << R"( <DataArray type="Float64" Name=")" << var << "\">"
799 << endl;
800 outfile << " ";
801 for (int n = 0; n < m_lines.size(); ++n)
802 {
803 const Array<OneD, NekDouble> phys =
804 m_phys + m_phys_offset[expansion] + n * npoints_per_line;
805
806 for (i = 0; i < nq; ++i)
807 {
808 outfile << (fabs(phys[i]) < NekConstants::kNekZeroTol ? 0 : phys[i])
809 << " ";
810 }
811 }
812 outfile << endl;
813 outfile << " </DataArray>" << endl;
814}
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
static const NekDouble kNekZeroTol

References Nektar::NekConstants::kNekZeroTol, m_lines, Nektar::MultiRegions::ExpList::m_phys, and Nektar::MultiRegions::ExpList::m_phys_offset.

Member Data Documentation

◆ m_dealiasing

bool Nektar::MultiRegions::ExpListHomogeneous2D::m_dealiasing
private

Definition at line 226 of file ExpListHomogeneous2D.h.

Referenced by ExpListHomogeneous2D().

◆ m_FFT_y

LibUtilities::NektarFFTSharedPtr Nektar::MultiRegions::ExpListHomogeneous2D::m_FFT_y

Definition at line 116 of file ExpListHomogeneous2D.h.

Referenced by ExpListHomogeneous2D(), and Homogeneous2DTrans().

◆ m_FFT_z

LibUtilities::NektarFFTSharedPtr Nektar::MultiRegions::ExpListHomogeneous2D::m_FFT_z

Definition at line 117 of file ExpListHomogeneous2D.h.

Referenced by ExpListHomogeneous2D(), and Homogeneous2DTrans().

◆ m_homogeneous2DBlockMat

Homo2DBlockMatrixMapShPtr Nektar::MultiRegions::ExpListHomogeneous2D::m_homogeneous2DBlockMat
protected

Definition at line 139 of file ExpListHomogeneous2D.h.

Referenced by GetHomogeneous2DBlockMatrix().

◆ m_homogeneousBasis_y

LibUtilities::BasisSharedPtr Nektar::MultiRegions::ExpListHomogeneous2D::m_homogeneousBasis_y
protected

◆ m_homogeneousBasis_z

LibUtilities::BasisSharedPtr Nektar::MultiRegions::ExpListHomogeneous2D::m_homogeneousBasis_z
protected

◆ m_lhom_y

NekDouble Nektar::MultiRegions::ExpListHomogeneous2D::m_lhom_y
protected

◆ m_lhom_z

NekDouble Nektar::MultiRegions::ExpListHomogeneous2D::m_lhom_z
protected

◆ m_lines

Array<OneD, ExpListSharedPtr> Nektar::MultiRegions::ExpListHomogeneous2D::m_lines
protected

Vector of ExpList, will be filled with ExpList1D.

Definition at line 141 of file ExpListHomogeneous2D.h.

Referenced by Nektar::MultiRegions::ContField3DHomogeneous2D::ContField3DHomogeneous2D(), Nektar::MultiRegions::DisContField3DHomogeneous2D::DisContField3DHomogeneous2D(), Nektar::MultiRegions::ExpList2DHomogeneous2D::ExpList2DHomogeneous2D(), Nektar::MultiRegions::ExpList3DHomogeneous2D::ExpList3DHomogeneous2D(), ExpListHomogeneous2D(), GenHomogeneous2DBlockMatrix(), Homogeneous2DTrans(), Nektar::MultiRegions::ExpList2DHomogeneous2D::SetCoeffPhys(), Nektar::MultiRegions::ExpList3DHomogeneous2D::SetCoeffPhys(), Nektar::MultiRegions::DisContField3DHomogeneous2D::SetupBoundaryConditions(), v_AppendFieldData(), v_BwdTrans(), Nektar::MultiRegions::ContField3DHomogeneous2D::v_ClearGlobalLinSysManager(), v_DealiasedProd(), Nektar::MultiRegions::DisContField3DHomogeneous2D::v_EvaluateBoundaryConditions(), v_ExtractDataToCoeffs(), v_FwdTrans(), v_FwdTransBndConstrained(), v_FwdTransLocalElmt(), Nektar::MultiRegions::DisContField3DHomogeneous2D::v_GetBoundaryToElmtMap(), Nektar::MultiRegions::ExpList3DHomogeneous2D::v_GetCoords(), Nektar::MultiRegions::ExpList2DHomogeneous2D::v_GetCoords(), v_GetFieldDefinitions(), v_GetNumElmts(), Nektar::MultiRegions::ContField3DHomogeneous2D::v_GlobalToLocal(), Nektar::MultiRegions::ContField3DHomogeneous2D::v_HelmSolve(), Nektar::MultiRegions::DisContField3DHomogeneous2D::v_HelmSolve(), Nektar::MultiRegions::ContField3DHomogeneous2D::v_ImposeDirichletConditions(), v_IProductWRTBase(), Nektar::MultiRegions::ExpList3DHomogeneous2D::v_L2(), Nektar::MultiRegions::ContField3DHomogeneous2D::v_LocalToGlobal(), v_PhysDeriv(), and v_WriteVtkPieceData().

◆ m_ny

int Nektar::MultiRegions::ExpListHomogeneous2D::m_ny
protected

◆ m_nz

int Nektar::MultiRegions::ExpListHomogeneous2D::m_nz
protected

◆ m_paddingBasis_y

LibUtilities::BasisSharedPtr Nektar::MultiRegions::ExpListHomogeneous2D::m_paddingBasis_y
protected

Base expansion in y direction.

Definition at line 134 of file ExpListHomogeneous2D.h.

Referenced by SetPaddingBase().

◆ m_paddingBasis_z

LibUtilities::BasisSharedPtr Nektar::MultiRegions::ExpListHomogeneous2D::m_paddingBasis_z
protected

Base expansion in z direction.

Definition at line 136 of file ExpListHomogeneous2D.h.

Referenced by SetPaddingBase().

◆ m_padsize_y

int Nektar::MultiRegions::ExpListHomogeneous2D::m_padsize_y
private

Definition at line 227 of file ExpListHomogeneous2D.h.

Referenced by SetPaddingBase(), and v_DealiasedProd().

◆ m_padsize_z

int Nektar::MultiRegions::ExpListHomogeneous2D::m_padsize_z
private

Definition at line 228 of file ExpListHomogeneous2D.h.

Referenced by SetPaddingBase(), and v_DealiasedProd().

◆ m_tmpIN

Array<OneD, NekDouble> Nektar::MultiRegions::ExpListHomogeneous2D::m_tmpIN

Definition at line 118 of file ExpListHomogeneous2D.h.

Referenced by Homogeneous2DTrans().

◆ m_tmpOUT

Array<OneD, NekDouble> Nektar::MultiRegions::ExpListHomogeneous2D::m_tmpOUT

Definition at line 119 of file ExpListHomogeneous2D.h.

Referenced by Homogeneous2DTrans().

◆ m_transposition

LibUtilities::TranspositionSharedPtr Nektar::MultiRegions::ExpListHomogeneous2D::m_transposition

◆ m_useFFT

bool Nektar::MultiRegions::ExpListHomogeneous2D::m_useFFT

◆ m_Ycomm

LibUtilities::CommSharedPtr Nektar::MultiRegions::ExpListHomogeneous2D::m_Ycomm

Definition at line 122 of file ExpListHomogeneous2D.h.

Referenced by ExpListHomogeneous2D().

◆ m_Zcomm

LibUtilities::CommSharedPtr Nektar::MultiRegions::ExpListHomogeneous2D::m_Zcomm

Definition at line 123 of file ExpListHomogeneous2D.h.

Referenced by ExpListHomogeneous2D().

◆ MatBwdPAD

DNekMatSharedPtr Nektar::MultiRegions::ExpListHomogeneous2D::MatBwdPAD
private

Definition at line 230 of file ExpListHomogeneous2D.h.

Referenced by SetPaddingBase().

◆ MatFwdPAD

DNekMatSharedPtr Nektar::MultiRegions::ExpListHomogeneous2D::MatFwdPAD
private

Definition at line 229 of file ExpListHomogeneous2D.h.

Referenced by SetPaddingBase().