Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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:
Inheritance graph
[legend]
Collaboration diagram for Nektar::MultiRegions::ExpListHomogeneous2D:
Collaboration graph
[legend]

Public Member Functions

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

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
 
- Public Attributes inherited from Nektar::MultiRegions::ExpList
ExpansionType m_expType
 

Protected Member Functions

DNekBlkMatSharedPtr GenHomogeneous2DBlockMatrix (Homogeneous2DMatType mattype, CoeffState coeffstate=eLocal) const
 
DNekBlkMatSharedPtr GetHomogeneous2DBlockMatrix (Homogeneous2DMatType mattype, CoeffState coeffstate=eLocal) const
 
virtual int v_GetNumElmts (void)
 
virtual void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
 
virtual void v_FwdTrans_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
 
virtual void v_BwdTrans_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
 
virtual void v_IProductWRTBase_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual std::vector
< LibUtilities::FieldDefinitionsSharedPtr
v_GetFieldDefinitions (void)
 
virtual void v_GetFieldDefinitions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef)
 
virtual void v_AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
 
virtual void v_AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, Array< OneD, NekDouble > &coeffs)
 
virtual void v_ExtractDataToCoeffs (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
 Extract data from raw field data into expansion list. More...
 
virtual void v_WriteVtkPieceData (std::ostream &outfile, int expansion, std::string var)
 
virtual void v_HomogeneousFwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
 
virtual void v_HomogeneousBwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
 
virtual void v_DealiasedProd (const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
 
virtual void v_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 (Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
 
- Protected Member Functions inherited from Nektar::MultiRegions::ExpList
boost::shared_ptr< DNekMatGenGlobalMatrixFull (const GlobalLinSysKey &mkey, const boost::shared_ptr< AssemblyMapCG > &locToGloMap)
 
const DNekScalBlkMatSharedPtr GenBlockMatrix (const GlobalMatrixKey &gkey)
 This function assembles the block diagonal matrix of local matrices of the type mtype. More...
 
const DNekScalBlkMatSharedPtrGetBlockMatrix (const GlobalMatrixKey &gkey)
 
void MultiplyByBlockMatrix (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
boost::shared_ptr< GlobalMatrixGenGlobalMatrix (const GlobalMatrixKey &mkey, const boost::shared_ptr< AssemblyMapCG > &locToGloMap)
 Generates a global matrix from the given key and map. More...
 
void GlobalEigenSystem (const boost::shared_ptr< DNekMat > &Gmat, Array< OneD, NekDouble > &EigValsReal, Array< OneD, NekDouble > &EigValsImag, Array< OneD, NekDouble > &EigVecs=NullNekDouble1DArray)
 
boost::shared_ptr< GlobalLinSysGenGlobalLinSys (const GlobalLinSysKey &mkey, const boost::shared_ptr< AssemblyMapCG > &locToGloMap)
 This operation constructs the global linear system of type mkey. More...
 
boost::shared_ptr< GlobalLinSysGenGlobalBndLinSys (const GlobalLinSysKey &mkey, const AssemblyMapSharedPtr &locToGloMap)
 Generate a GlobalLinSys from information provided by the key "mkey" and the mapping provided in LocToGloBaseMap. More...
 
void ReadGlobalOptimizationParameters ()
 
virtual const Array< OneD,
const boost::shared_ptr
< ExpList > > & 
v_GetBndCondExpansions (void)
 
virtual boost::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 boost::shared_ptr
< ExpList > & 
v_GetTrace ()
 
virtual boost::shared_ptr
< AssemblyMapDG > & 
v_GetTraceMap ()
 
virtual const Array< OneD,
const int > & 
v_GetTraceBndMap ()
 
virtual void v_GetNormals (Array< OneD, Array< OneD, NekDouble > > &normals)
 
virtual void v_AddTraceIntegral (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddTraceIntegral (const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 
virtual void v_AddFwdBwdTraceIntegral (const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
 
virtual void v_GetFwdBwdTracePhys (Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
 
virtual void v_GetFwdBwdTracePhys (const Array< OneD, const NekDouble > &field, Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
 
virtual const 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, CoeffState coeffstate)
 
virtual void v_HelmSolve (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const FlagList &flags, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff, const Array< OneD, const NekDouble > &dirForcing)
 
virtual void v_LinearAdvectionDiffusionReactionSolve (const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 
virtual void v_LinearAdvectionReactionSolve (const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 
virtual void v_ImposeDirichletConditions (Array< OneD, NekDouble > &outarray)
 
virtual void v_FillBndCondFromField ()
 
virtual void v_Reset ()
 Reset geometry information, metrics, matrix managers and geometry information. More...
 
virtual void v_LocalToGlobal (void)
 
virtual void v_GlobalToLocal (void)
 
virtual void v_SmoothField (Array< OneD, NekDouble > &field)
 
virtual void v_GeneralMatrixOp (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
 
virtual void v_GetCoords (Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
 
virtual void v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
 
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 ()
 
virtual void v_GetBoundaryToElmtMap (Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
 
virtual void v_GetBndElmtExpansion (int i, boost::shared_ptr< ExpList > &result)
 
virtual void v_ExtractElmtToBndPhys (int i, Array< OneD, NekDouble > &elmt, Array< OneD, NekDouble > &boundary)
 
virtual void v_ExtractPhysToBndElmt (int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bndElmt)
 
virtual void v_GetBoundaryNormals (int i, Array< OneD, Array< OneD, NekDouble > > &normals)
 
virtual void v_ReadGlobalOptimizationParameters ()
 
virtual void v_ExtractCoeffsToCoeffs (const boost::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
 
virtual void v_WriteTecplotHeader (std::ostream &outfile, std::string var="")
 
virtual void v_WriteTecplotZone (std::ostream &outfile, int expansion)
 
virtual void v_WriteTecplotField (std::ostream &outfile, int expansion)
 
virtual void v_WriteTecplotConnectivity (std::ostream &outfile, int expansion)
 
virtual void v_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 Array< OneD, const
NekDouble
v_HomogeneousEnergy (void)
 
virtual
LibUtilities::TranspositionSharedPtr 
v_GetTransposition (void)
 
virtual NekDouble v_GetHomoLen (void)
 
virtual Array< OneD, const
unsigned int > 
v_GetZIDs (void)
 
virtual Array< OneD, const
unsigned int > 
v_GetYIDs (void)
 
virtual void v_PhysInterp1DScaled (const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_PhysGalerkinProjection1DScaled (const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
virtual void v_ClearGlobalLinSysManager (void)
 
void ExtractFileBCs (const std::string &fileName, const std::string &varName, const boost::shared_ptr< ExpList > locExpList)
 

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
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
LibUtilities::SessionReaderSharedPtr m_session
 Session. More...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 Mesh associated with this expansion list. More...
 
int m_ncoeffs
 The total number of local degrees of freedom. m_ncoeffs $=N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_l$. More...
 
int m_npoints
 
Array< OneD, NekDoublem_coeffs
 Concatenation of all local expansion coefficients. More...
 
Array< OneD, NekDoublem_phys
 The global expansion evaluated at the quadrature points. More...
 
bool m_physState
 The state of the array m_phys. More...
 
boost::shared_ptr
< LocalRegions::ExpansionVector
m_exp
 The list of local expansions. More...
 
Collections::CollectionVector m_collections
 
std::vector< int > m_coll_coeff_offset
 Offset of elemental data into the array m_coeffs. More...
 
std::vector< int > m_coll_phys_offset
 Offset of elemental data into the array m_phys. More...
 
Array< OneD, int > m_coeff_offset
 Offset of elemental data into the array m_coeffs. More...
 
Array< OneD, int > m_phys_offset
 Offset of elemental data into the array m_phys. More...
 
Array< OneD, int > m_offset_elmt_id
 Array containing the element id m_offset_elmt_id[n] that the n^th consecutive block of data in m_coeffs and m_phys is associated, i.e. for an array of constant expansion size and single shape elements m_phys[n*m_npoints] is the data related to m_exp[m_offset_elmt_id[n]];. More...
 
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
 
BlockMatrixMapShPtr m_blockMat
 
bool m_WaveSpace
 
boost::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 81 of file ExpListHomogeneous2D.h.

Constructor & Destructor Documentation

Nektar::MultiRegions::ExpListHomogeneous2D::ExpListHomogeneous2D ( )

Default constructor.

Definition at line 49 of file ExpListHomogeneous2D.cpp.

49  :
50  ExpList(),
53  m_lhom_y(1),
54  m_lhom_z(1),
55  m_homogeneous2DBlockMat(MemoryManager<Homo2DBlockMatrixMap>::AllocateSharedPtr())
56  {
57  }
ExpList()
The default constructor.
Definition: ExpList.cpp:95
static BasisSharedPtr NullBasisSharedPtr
Definition: Basis.h:358
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.
LibUtilities::BasisSharedPtr m_homogeneousBasis_y
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
Nektar::MultiRegions::ExpListHomogeneous2D::ExpListHomogeneous2D ( 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 59 of file ExpListHomogeneous2D.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, ASSERTL2, Nektar::LibUtilities::BasisManager(), Nektar::LibUtilities::NekFactory< tKey, tBase, >::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, m_transposition, m_useFFT, m_Ycomm, m_Zcomm, Nektar::LibUtilities::NullBasisKey(), and SetPaddingBase().

65  :
66  ExpList(pSession),
67  m_useFFT(useFFT),
68  m_lhom_y(lhom_y),
69  m_lhom_z(lhom_z),
70  m_homogeneous2DBlockMat(MemoryManager<Homo2DBlockMatrixMap>::AllocateSharedPtr()),
71  m_dealiasing(dealiasing)
72  {
73  ASSERTL2(HomoBasis_y != LibUtilities::NullBasisKey,
74  "Homogeneous Basis in y direction is a null basis");
75  ASSERTL2(HomoBasis_z != LibUtilities::NullBasisKey,
76  "Homogeneous Basis in z direction is a null basis");
77 
80 
82 
83  m_Ycomm = m_comm->GetColumnComm()->GetRowComm();
84  m_Zcomm = m_comm->GetColumnComm()->GetRowComm();
85 
86  m_ny = m_homogeneousBasis_y->GetNumPoints()/m_Ycomm->GetSize();
87  m_nz = m_homogeneousBasis_z->GetNumPoints()/m_Zcomm->GetSize();
88 
89  m_lines = Array<OneD,ExpListSharedPtr>(m_ny*m_nz);
90 
91  if(m_useFFT)
92  {
95  }
96 
97  if(m_dealiasing)
98  {
99  ASSERTL0(m_comm->GetColumnComm()->GetSize() == 1,"Remove dealiasing if you want to run in parallel");
100  SetPaddingBase();
101  }
102  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
LibUtilities::TranspositionSharedPtr m_transposition
ExpList()
The default constructor.
Definition: ExpList.cpp:95
tBaseSharedPtr CreateInstance(tKey idKey BOOST_PP_COMMA_IF(MAX_PARAM) BOOST_PP_ENUM_BINARY_PARAMS(MAX_PARAM, tParam, x))
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:162
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
LibUtilities::NektarFFTSharedPtr m_FFT_z
NekDouble m_lhom_z
Width of homogeneous direction z.
NektarFFTFactory & GetNektarFFTFactory()
Definition: NektarFFT.cpp:69
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
BasisManagerT & BasisManager(void)
LibUtilities::BasisSharedPtr m_homogeneousBasis_z
Base expansion in z direction.
NekDouble m_lhom_y
Width of homogeneous direction y.
int m_ny
Number of modes = number of poitns in y direction.
LibUtilities::BasisSharedPtr m_homogeneousBasis_y
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
int m_nz
Number of modes = number of poitns in z direction.
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: ExpList.h:907
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:240
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.
LibUtilities::NektarFFTSharedPtr m_FFT_y
Nektar::MultiRegions::ExpListHomogeneous2D::ExpListHomogeneous2D ( const ExpListHomogeneous2D In)

Copy constructor.

Parameters
InExpListHomogeneous2D object to copy.

Definition at line 108 of file ExpListHomogeneous2D.cpp.

References m_lines.

108  :
109  ExpList(In,false),
110  m_useFFT(In.m_useFFT),
111  m_FFT_y(In.m_FFT_y),
112  m_FFT_z(In.m_FFT_z),
113  m_transposition(In.m_transposition),
114  m_Ycomm(In.m_Ycomm),
115  m_Zcomm(In.m_Ycomm),
116  m_homogeneousBasis_y(In.m_homogeneousBasis_y),
117  m_homogeneousBasis_z(In.m_homogeneousBasis_z),
118  m_lhom_y(In.m_lhom_y),
119  m_lhom_z(In.m_lhom_z),
120  m_homogeneous2DBlockMat(In.m_homogeneous2DBlockMat),
121  m_ny(In.m_ny),
122  m_nz(In.m_nz),
123  m_dealiasing(In.m_dealiasing),
124  m_padsize_y(In.m_padsize_y),
125  m_padsize_z(In.m_padsize_z),
126  MatFwdPAD(In.MatFwdPAD),
127  MatBwdPAD(In.MatBwdPAD)
128  {
129  m_lines = Array<OneD, ExpListSharedPtr>(In.m_lines.num_elements());
130  }
LibUtilities::TranspositionSharedPtr m_transposition
ExpList()
The default constructor.
Definition: ExpList.cpp:95
LibUtilities::NektarFFTSharedPtr m_FFT_z
NekDouble m_lhom_z
Width of homogeneous direction z.
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
LibUtilities::BasisSharedPtr m_homogeneousBasis_z
Base expansion in z direction.
NekDouble m_lhom_y
Width of homogeneous direction y.
int m_ny
Number of modes = number of poitns in y direction.
LibUtilities::BasisSharedPtr m_homogeneousBasis_y
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
int m_nz
Number of modes = number of poitns in z direction.
LibUtilities::NektarFFTSharedPtr m_FFT_y
Nektar::MultiRegions::ExpListHomogeneous2D::~ExpListHomogeneous2D ( )
virtual

Destructor.

Destructor

Definition at line 135 of file ExpListHomogeneous2D.cpp.

136  {
137  }

Member Function Documentation

void Nektar::MultiRegions::ExpListHomogeneous2D::DealiasedProd ( const Array< OneD, NekDouble > &  inarray1,
const Array< OneD, NekDouble > &  inarray2,
Array< OneD, NekDouble > &  outarray,
CoeffState  coeffstate = eLocal 
)
inline

Definition at line 257 of file ExpListHomogeneous2D.h.

References v_DealiasedProd().

261  {
262  v_DealiasedProd(inarray1,inarray2,outarray,coeffstate);
263  }
virtual void v_DealiasedProd(const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
DNekBlkMatSharedPtr Nektar::MultiRegions::ExpListHomogeneous2D::GenHomogeneous2DBlockMatrix ( Homogeneous2DMatType  mattype,
CoeffState  coeffstate = eLocal 
) const
protected

Definition at line 495 of file ExpListHomogeneous2D.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), 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, m_homogeneousBasis_y, m_homogeneousBasis_z, and m_lines.

Referenced by GetHomogeneous2DBlockMatrix().

496  {
497  int i;
498  int n_exp = 0;
499 
500  DNekMatSharedPtr loc_mat;
501  DNekBlkMatSharedPtr BlkMatrix;
502 
504 
505  int NumPoints = 0;
506  int NumModes = 0;
507  int NumPencils = 0;
508 
509  if((mattype == eForwardsCoeffSpaceY1D) || (mattype == eBackwardsCoeffSpaceY1D)
510  ||(mattype == eForwardsPhysSpaceY1D) || (mattype == eBackwardsPhysSpaceY1D))
511  {
512  Basis = m_homogeneousBasis_y;
513  NumPoints = m_homogeneousBasis_y->GetNumModes();
514  NumModes = m_homogeneousBasis_y->GetNumPoints();
515  NumPencils = m_homogeneousBasis_z->GetNumPoints();
516  }
517  else
518  {
519  Basis = m_homogeneousBasis_z;
520  NumPoints = m_homogeneousBasis_z->GetNumModes();
521  NumModes = m_homogeneousBasis_z->GetNumPoints();
522  NumPencils = m_homogeneousBasis_y->GetNumPoints();
523  }
524 
525  if((mattype == eForwardsCoeffSpaceY1D) || (mattype == eForwardsCoeffSpaceZ1D)
526  ||(mattype == eBackwardsCoeffSpaceY1D)||(mattype == eBackwardsCoeffSpaceZ1D))
527  {
528  n_exp = m_lines[0]->GetNcoeffs();
529  }
530  else
531  {
532  n_exp = m_lines[0]->GetTotPoints(); // will operatore on m_phys
533  }
534 
535  Array<OneD,unsigned int> nrows(n_exp);
536  Array<OneD,unsigned int> ncols(n_exp);
537 
538  if((mattype == eForwardsCoeffSpaceY1D)||(mattype == eForwardsPhysSpaceY1D) ||
539  (mattype == eForwardsCoeffSpaceZ1D)||(mattype == eForwardsPhysSpaceZ1D))
540  {
541  nrows = Array<OneD, unsigned int>(n_exp*NumPencils,NumModes);
542  ncols = Array<OneD, unsigned int>(n_exp*NumPencils,NumPoints);
543  }
544  else
545  {
546  nrows = Array<OneD, unsigned int>(n_exp*NumPencils,NumPoints);
547  ncols = Array<OneD, unsigned int>(n_exp*NumPencils,NumModes);
548  }
549 
550  MatrixStorage blkmatStorage = eDIAGONAL;
551  BlkMatrix = MemoryManager<DNekBlkMat>::AllocateSharedPtr(nrows,ncols,blkmatStorage);
552 
553  StdRegions::StdSegExp StdSeg(Basis->GetBasisKey());
554 
555  if((mattype == eForwardsCoeffSpaceY1D)||(mattype == eForwardsPhysSpaceY1D) ||
556  (mattype == eForwardsCoeffSpaceZ1D)||(mattype == eForwardsPhysSpaceZ1D))
557  {
558  StdRegions::StdMatrixKey matkey(StdRegions::eFwdTrans,
559  StdSeg.DetShapeType(),
560  StdSeg);
561 
562  loc_mat = StdSeg.GetStdMatrix(matkey);
563  }
564  else
565  {
566  StdRegions::StdMatrixKey matkey(StdRegions::eBwdTrans,
567  StdSeg.DetShapeType(),
568  StdSeg);
569 
570  loc_mat = StdSeg.GetStdMatrix(matkey);
571  }
572 
573  // set up array of block matrices.
574  for(i = 0; i < (n_exp*NumPencils); ++i)
575  {
576  BlkMatrix->SetBlock(i,i,loc_mat);
577  }
578 
579  return BlkMatrix;
580  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
LibUtilities::BasisSharedPtr m_homogeneousBasis_z
Base expansion in z direction.
LibUtilities::BasisSharedPtr m_homogeneousBasis_y
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
boost::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:72
boost::shared_ptr< Basis > BasisSharedPtr
DNekBlkMatSharedPtr Nektar::MultiRegions::ExpListHomogeneous2D::GetHomogeneous2DBlockMatrix ( Homogeneous2DMatType  mattype,
CoeffState  coeffstate = eLocal 
) const
protected

Definition at line 480 of file ExpListHomogeneous2D.cpp.

References GenHomogeneous2DBlockMatrix(), Nektar::iterator, and m_homogeneous2DBlockMat.

Referenced by Homogeneous2DTrans().

481  {
482  Homo2DBlockMatrixMap::iterator matrixIter = m_homogeneous2DBlockMat->find(mattype);
483 
484  if(matrixIter == m_homogeneous2DBlockMat->end())
485  {
486  return ((*m_homogeneous2DBlockMat)[mattype] =
487  GenHomogeneous2DBlockMatrix(mattype,coeffstate));
488  }
489  else
490  {
491  return matrixIter->second;
492  }
493  }
DNekBlkMatSharedPtr GenHomogeneous2DBlockMatrix(Homogeneous2DMatType mattype, CoeffState coeffstate=eLocal) const
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
void Nektar::MultiRegions::ExpListHomogeneous2D::Homogeneous2DTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
bool  IsForwards,
CoeffState  coeffstate = eLocal,
bool  Shuff = true,
bool  UnShuff = true 
)

Definition at line 357 of file ExpListHomogeneous2D.cpp.

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, and m_useFFT.

Referenced by v_HomogeneousBwdTrans(), and v_HomogeneousFwdTrans().

363  {
364  if(m_useFFT)
365  {
366 
367  int n = m_lines.num_elements(); //number of Fourier points in the Fourier directions (x-z grid)
368  int s = inarray.num_elements(); //number of total points = n. of Fourier points * n. of points per line
369  int p = s/n; //number of points per line = n of Fourier transform required
370 
371  Array<OneD, NekDouble> fft_in(s);
372  Array<OneD, NekDouble> fft_out(s);
373 
374  m_transposition->Transpose(inarray,fft_in,false,LibUtilities::eXtoYZ);
375 
376  if(IsForwards)
377  {
378  for(int i=0;i<(p*m_nz);i++)
379  {
380  m_FFT_y->FFTFwdTrans(m_tmpIN = fft_in + i*m_ny, m_tmpOUT = fft_out + i*m_ny);
381  }
382 
383  }
384  else
385  {
386  for(int i=0;i<(p*m_nz);i++)
387  {
388  m_FFT_y->FFTBwdTrans(m_tmpIN = fft_in + i*m_ny, m_tmpOUT = fft_out + i*m_ny);
389  }
390  }
391 
392  m_transposition->Transpose(fft_out,fft_in,false,LibUtilities::eYZtoZY);
393 
394  if(IsForwards)
395  {
396  for(int i=0;i<(p*m_ny);i++)
397  {
398  m_FFT_z->FFTFwdTrans(m_tmpIN = fft_in + i*m_nz, m_tmpOUT = fft_out + i*m_nz);
399  }
400 
401  }
402  else
403  {
404  for(int i=0;i<(p*m_ny);i++)
405  {
406  m_FFT_z->FFTBwdTrans(m_tmpIN = fft_in + i*m_nz, m_tmpOUT = fft_out + i*m_nz);
407  }
408  }
409 
410  //TODO: required ZYtoX routine
411  m_transposition->Transpose(fft_out,fft_in,false,LibUtilities::eZYtoYZ);
412 
413  m_transposition->Transpose(fft_in,outarray,false,LibUtilities::eYZtoX);
414 
415  }
416  else
417  {
418  DNekBlkMatSharedPtr blkmatY;
419  DNekBlkMatSharedPtr blkmatZ;
420 
421  if(inarray.num_elements() == m_npoints) //transform phys space
422  {
423  if(IsForwards)
424  {
427  }
428  else
429  {
432  }
433  }
434  else
435  {
436  if(IsForwards)
437  {
440  }
441  else
442  {
445  }
446  }
447 
448  int nrowsY = blkmatY->GetRows();
449  int ncolsY = blkmatY->GetColumns();
450 
451  Array<OneD, NekDouble> sortedinarrayY(ncolsY);
452  Array<OneD, NekDouble> sortedoutarrayY(nrowsY);
453 
454  int nrowsZ = blkmatZ->GetRows();
455  int ncolsZ = blkmatZ->GetColumns();
456 
457  Array<OneD, NekDouble> sortedinarrayZ(ncolsZ);
458  Array<OneD, NekDouble> sortedoutarrayZ(nrowsZ);
459 
460  NekVector<NekDouble> inY (ncolsY,sortedinarrayY,eWrapper);
461  NekVector<NekDouble> outY(nrowsY,sortedoutarrayY,eWrapper);
462 
463  NekVector<NekDouble> inZ (ncolsZ,sortedinarrayZ,eWrapper);
464  NekVector<NekDouble> outZ(nrowsZ,sortedoutarrayZ,eWrapper);
465 
466  m_transposition->Transpose(inarray,sortedinarrayY,!IsForwards,LibUtilities::eXtoYZ);
467 
468  outY = (*blkmatY)*inY;
469 
470  m_transposition->Transpose(sortedoutarrayY,sortedinarrayZ,false,LibUtilities::eYZtoZY);
471 
472  outZ = (*blkmatZ)*inZ;
473 
474  m_transposition->Transpose(sortedoutarrayZ,sortedoutarrayY,false,LibUtilities::eZYtoYZ);
475 
476  m_transposition->Transpose(sortedoutarrayY,outarray,false,LibUtilities::eYZtoX);
477  }
478  }
LibUtilities::TranspositionSharedPtr m_transposition
LibUtilities::NektarFFTSharedPtr m_FFT_z
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
int m_ny
Number of modes = number of poitns in y direction.
int m_nz
Number of modes = number of poitns in z direction.
boost::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:72
LibUtilities::NektarFFTSharedPtr m_FFT_y
DNekBlkMatSharedPtr GetHomogeneous2DBlockMatrix(Homogeneous2DMatType mattype, CoeffState coeffstate=eLocal) const
void Nektar::MultiRegions::ExpListHomogeneous2D::HomogeneousBwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
CoeffState  coeffstate = eLocal,
bool  Shuff = true,
bool  UnShuff = true 
)
inline

Definition at line 248 of file ExpListHomogeneous2D.h.

References v_HomogeneousBwdTrans().

Referenced by v_BwdTrans(), v_BwdTrans_IterPerExp(), v_DealiasedProd(), and v_PhysDeriv().

253  {
254  v_HomogeneousBwdTrans(inarray,outarray,coeffstate,Shuff,UnShuff);
255  }
virtual void v_HomogeneousBwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
void Nektar::MultiRegions::ExpListHomogeneous2D::HomogeneousFwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
CoeffState  coeffstate = eLocal,
bool  Shuff = true,
bool  UnShuff = true 
)
inline

Definition at line 239 of file ExpListHomogeneous2D.h.

References v_HomogeneousFwdTrans().

Referenced by v_DealiasedProd(), v_FwdTrans(), v_FwdTrans_IterPerExp(), Nektar::MultiRegions::ContField3DHomogeneous2D::v_HelmSolve(), Nektar::MultiRegions::DisContField3DHomogeneous2D::v_HelmSolve(), and v_PhysDeriv().

244  {
245  v_HomogeneousFwdTrans(inarray,outarray,coeffstate,Shuff,UnShuff);
246  }
virtual void v_HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
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 957 of file ExpListHomogeneous2D.cpp.

References v_PhysDeriv().

962  {
963  v_PhysDeriv(inarray,out_d0,out_d1,out_d2);
964  }
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)
void Nektar::MultiRegions::ExpListHomogeneous2D::PhysDeriv ( Direction  edir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d 
)

Definition at line 966 of file ExpListHomogeneous2D.cpp.

References v_PhysDeriv().

969  {
970  //convert int into enum
971  v_PhysDeriv(edir,inarray,out_d);
972  }
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)
void Nektar::MultiRegions::ExpListHomogeneous2D::SetPaddingBase ( void  )

Definition at line 974 of file ExpListHomogeneous2D.cpp.

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

Referenced by ExpListHomogeneous2D().

975  {
976  NekDouble size_y = 1.5*m_ny;
977  NekDouble size_z = 1.5*m_nz;
978  m_padsize_y = int(size_y);
979  m_padsize_z = int(size_z);
980 
981  const LibUtilities::PointsKey Ppad_y(m_padsize_y,LibUtilities::eFourierEvenlySpaced);
982  const LibUtilities::BasisKey Bpad_y(LibUtilities::eFourier,m_padsize_y,Ppad_y);
983 
984  const LibUtilities::PointsKey Ppad_z(m_padsize_z,LibUtilities::eFourierEvenlySpaced);
985  const LibUtilities::BasisKey Bpad_z(LibUtilities::eFourier,m_padsize_z,Ppad_z);
986 
989 
990  StdRegions::StdQuadExp StdQuad(m_paddingBasis_y->GetBasisKey(),m_paddingBasis_z->GetBasisKey());
991 
992  StdRegions::StdMatrixKey matkey1(StdRegions::eFwdTrans,StdQuad.DetShapeType(),StdQuad);
993  StdRegions::StdMatrixKey matkey2(StdRegions::eBwdTrans,StdQuad.DetShapeType(),StdQuad);
994 
995  MatFwdPAD = StdQuad.GetStdMatrix(matkey1);
996  MatBwdPAD = StdQuad.GetStdMatrix(matkey2);
997  }
LibUtilities::BasisSharedPtr m_paddingBasis_z
Base expansion in z direction.
Fourier Expansion .
Definition: BasisType.h:52
BasisManagerT & BasisManager(void)
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:64
double NekDouble
int m_ny
Number of modes = number of poitns in y direction.
int m_nz
Number of modes = number of poitns in z direction.
LibUtilities::BasisSharedPtr m_paddingBasis_y
Base expansion in y direction.
void Nektar::MultiRegions::ExpListHomogeneous2D::v_AppendFieldData ( LibUtilities::FieldDefinitionsSharedPtr fielddef,
std::vector< NekDouble > &  fielddata 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 681 of file ExpListHomogeneous2D.cpp.

References Nektar::MultiRegions::ExpList::m_coeffs.

682  {
683  v_AppendFieldData(fielddef,fielddata,m_coeffs);
684  }
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:939
virtual void v_AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
void Nektar::MultiRegions::ExpListHomogeneous2D::v_AppendFieldData ( LibUtilities::FieldDefinitionsSharedPtr fielddef,
std::vector< NekDouble > &  fielddata,
Array< OneD, NekDouble > &  coeffs 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 652 of file ExpListHomogeneous2D.cpp.

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

653  {
654  int i,k;
655 
656  int NumMod_y = m_homogeneousBasis_y->GetNumModes();
657  int NumMod_z = m_homogeneousBasis_z->GetNumModes();
658 
659  int ncoeffs_per_line = m_lines[0]->GetNcoeffs();
660 
661  // Determine mapping from element ids to location in
662  // expansion list
663  map<int, int> ElmtID_to_ExpID;
664  for(i = 0; i < m_lines[0]->GetExpSize(); ++i)
665  {
666  ElmtID_to_ExpID[(*m_exp)[i]->GetGeom()->GetGlobalID()] = i;
667  }
668 
669  for(i = 0; i < fielddef->m_elementIDs.size(); ++i)
670  {
671  int eid = ElmtID_to_ExpID[fielddef->m_elementIDs[i]];
672  int datalen = (*m_exp)[eid]->GetNcoeffs();
673 
674  for(k = 0; k < (NumMod_y*NumMod_z); ++k)
675  {
676  fielddata.insert(fielddata.end(),&coeffs[m_coeff_offset[eid]+k*ncoeffs_per_line],&coeffs[m_coeff_offset[eid]+k*ncoeffs_per_line]+datalen);
677  }
678  }
679  }
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:988
LibUtilities::BasisSharedPtr m_homogeneousBasis_z
Base expansion in z direction.
LibUtilities::BasisSharedPtr m_homogeneousBasis_y
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
void Nektar::MultiRegions::ExpListHomogeneous2D::v_BwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
CoeffState  coeffstate 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 288 of file ExpListHomogeneous2D.cpp.

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

289  {
290  int cnt = 0, cnt1 = 0;
291  Array<OneD, NekDouble> tmparray;
292  int nlines = m_lines.num_elements();
293 
294  for(int n = 0; n < nlines; ++n)
295  {
296  m_lines[n]->BwdTrans(inarray+cnt, tmparray = outarray + cnt1,
297  coeffstate);
298  cnt += m_lines[n]->GetNcoeffs();
299  cnt1 += m_lines[n]->GetTotPoints();
300  }
301  if(!m_WaveSpace)
302  {
303  HomogeneousBwdTrans(outarray,outarray);
304  }
305  }
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
void HomogeneousBwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
void Nektar::MultiRegions::ExpListHomogeneous2D::v_BwdTrans_IterPerExp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

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 307 of file ExpListHomogeneous2D.cpp.

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

308  {
309  int cnt = 0, cnt1 = 0;
310  Array<OneD, NekDouble> tmparray;
311  int nlines = m_lines.num_elements();
312 
313  for(int n = 0; n < nlines; ++n)
314  {
315  m_lines[n]->BwdTrans_IterPerExp(inarray+cnt, tmparray = outarray + cnt1);
316 
317  cnt += m_lines[n]->GetNcoeffs();
318  cnt1 += m_lines[n]->GetTotPoints();
319  }
320  if(!m_WaveSpace)
321  {
322  HomogeneousBwdTrans(outarray,outarray);
323  }
324  }
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
void HomogeneousBwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
void Nektar::MultiRegions::ExpListHomogeneous2D::v_DealiasedProd ( const Array< OneD, NekDouble > &  inarray1,
const Array< OneD, NekDouble > &  inarray2,
Array< OneD, NekDouble > &  outarray,
CoeffState  coeffstate = eLocal 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 159 of file ExpListHomogeneous2D.cpp.

References Nektar::eWrapper, Nektar::LibUtilities::eXtoYZ, Nektar::LibUtilities::eYZtoX, HomogeneousBwdTrans(), 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().

Referenced by DealiasedProd().

163  {
164  int npoints = outarray.num_elements(); // number of total physical points
165  int nlines = m_lines.num_elements(); // number of lines == number of Fourier modes = number of Fourier coeff = number of points per slab
166  int nslabs = npoints/nlines; // number of slabs = numebr of physical points per line
167 
168  Array<OneD, NekDouble> V1(npoints);
169  Array<OneD, NekDouble> V2(npoints);
170  Array<OneD, NekDouble> V1V2(npoints);
171  Array<OneD, NekDouble> ShufV1(npoints);
172  Array<OneD, NekDouble> ShufV2(npoints);
173  Array<OneD, NekDouble> ShufV1V2(npoints);
174 
175  if(m_WaveSpace)
176  {
177  V1 = inarray1;
178  V2 = inarray2;
179  }
180  else
181  {
182  HomogeneousFwdTrans(inarray1,V1,coeffstate);
183  HomogeneousFwdTrans(inarray2,V2,coeffstate);
184  }
185 
186  m_transposition->Transpose(V1,ShufV1,false,LibUtilities::eXtoYZ);
187  m_transposition->Transpose(V2,ShufV2,false,LibUtilities::eXtoYZ);
188 
189  Array<OneD, NekDouble> PadV1_slab_coeff(m_padsize_y*m_padsize_z,0.0);
190  Array<OneD, NekDouble> PadV2_slab_coeff(m_padsize_y*m_padsize_z,0.0);
191  Array<OneD, NekDouble> PadRe_slab_coeff(m_padsize_y*m_padsize_z,0.0);
192 
193  Array<OneD, NekDouble> PadV1_slab_phys(m_padsize_y*m_padsize_z,0.0);
194  Array<OneD, NekDouble> PadV2_slab_phys(m_padsize_y*m_padsize_z,0.0);
195  Array<OneD, NekDouble> PadRe_slab_phys(m_padsize_y*m_padsize_z,0.0);
196 
197  NekVector<NekDouble> PadIN_V1(m_padsize_y*m_padsize_z,PadV1_slab_coeff,eWrapper);
198  NekVector<NekDouble> PadOUT_V1(m_padsize_y*m_padsize_z,PadV1_slab_phys,eWrapper);
199 
200  NekVector<NekDouble> PadIN_V2(m_padsize_y*m_padsize_z,PadV2_slab_coeff,eWrapper);
201  NekVector<NekDouble> PadOUT_V2(m_padsize_y*m_padsize_z,PadV2_slab_phys,eWrapper);
202 
203  NekVector<NekDouble> PadIN_Re(m_padsize_y*m_padsize_z,PadRe_slab_phys,eWrapper);
204  NekVector<NekDouble> PadOUT_Re(m_padsize_y*m_padsize_z,PadRe_slab_coeff,eWrapper);
205 
206  //Looping on the slabs
207  for(int j = 0 ; j< nslabs ; j++)
208  {
209  //Copying the j-th slab of size N*M into a bigger slab of lenght 2*N*M
210  //We are in Fourier space
211  for(int i = 0 ; i< m_nz ; i++)
212  {
213  Vmath::Vcopy(m_ny,&(ShufV1[i*m_ny + j*nlines]),1,&(PadV1_slab_coeff[i*2*m_ny]),1);
214  Vmath::Vcopy(m_ny,&(ShufV2[i*m_ny + j*nlines]),1,&(PadV2_slab_coeff[i*2*m_ny]),1);
215  }
216 
217  //Moving to physical space using the padded system
218  PadOUT_V1 = (*MatBwdPAD)*PadIN_V1;
219  PadOUT_V2 = (*MatBwdPAD)*PadIN_V2;
220 
221  //Perfroming the vectors multiplication in physical
222  //space on the padded system
223  Vmath::Vmul(m_padsize_y*m_padsize_z,PadV1_slab_phys,1,PadV2_slab_phys,1,PadRe_slab_phys,1);
224 
225  //Moving back the result (V1*V2)_phys in Fourier
226  //space, padded system
227  PadOUT_Re = (*MatFwdPAD)*PadIN_Re;
228 
229  //Copying the first half of the padded pencil in the
230  //full vector (Fourier space)
231  for (int i = 0; i < m_nz; i++)
232  {
233  Vmath::Vcopy(m_ny,&(PadRe_slab_coeff[i*2*m_ny]),1,&(ShufV1V2[i*m_ny + j*nlines]),1);
234  }
235  }
236 
237  if(m_WaveSpace)
238  {
239  m_transposition->Transpose(ShufV1V2,outarray,false,LibUtilities::eYZtoX);
240  }
241  else
242  {
243  m_transposition->Transpose(ShufV1V2,V1V2,false,LibUtilities::eYZtoX);
244 
245  //Moving the results in physical space for the output
246  HomogeneousBwdTrans(V1V2,outarray,coeffstate);
247  }
248  }
LibUtilities::TranspositionSharedPtr m_transposition
void HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
int m_ny
Number of modes = number of poitns in y direction.
int m_nz
Number of modes = number of poitns in z direction.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
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.cpp:169
void HomogeneousBwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
void Nektar::MultiRegions::ExpListHomogeneous2D::v_ExtractDataToCoeffs ( LibUtilities::FieldDefinitionsSharedPtr fielddef,
std::vector< NekDouble > &  fielddata,
std::string &  field,
Array< OneD, NekDouble > &  coeffs 
)
protectedvirtual

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 687 of file ExpListHomogeneous2D.cpp.

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

688  {
689  int i,k;
690  int offset = 0;
691  int datalen = fielddata.size()/fielddef->m_fields.size();
692  int ncoeffs_per_line = m_lines[0]->GetNcoeffs();
693  int NumMod_y = m_homogeneousBasis_y->GetNumModes();
694  int NumMod_z = m_homogeneousBasis_z->GetNumModes();
695 
696  // Find data location according to field definition
697  for(i = 0; i < fielddef->m_fields.size(); ++i)
698  {
699  if(fielddef->m_fields[i] == field)
700  {
701  break;
702  }
703  offset += datalen;
704  }
705 
706  ASSERTL0(i!= fielddef->m_fields.size(),"Field not found in data file");
707 
708  // Determine mapping from element ids to location in
709  // expansion list
710  map<int, int> ElmtID_to_ExpID;
711  for(i = 0; i < m_lines[0]->GetExpSize(); ++i)
712  {
713  ElmtID_to_ExpID[(*m_exp)[i]->GetGeom()->GetGlobalID()] = i;
714  }
715 
716  for(i = 0; i < fielddef->m_elementIDs.size(); ++i)
717  {
718  int eid = ElmtID_to_ExpID[fielddef->m_elementIDs[i]];
719  int datalen = (*m_exp)[eid]->GetNcoeffs();
720 
721  for(k = 0; k < (NumMod_y*NumMod_z); ++k)
722  {
723  Vmath::Vcopy(datalen,&fielddata[offset],1,&coeffs[m_coeff_offset[eid] + k*ncoeffs_per_line],1);
724  offset += datalen;
725  }
726  }
727  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:988
LibUtilities::BasisSharedPtr m_homogeneousBasis_z
Base expansion in z direction.
LibUtilities::BasisSharedPtr m_homogeneousBasis_y
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
void Nektar::MultiRegions::ExpListHomogeneous2D::v_FwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
CoeffState  coeffstate 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 250 of file ExpListHomogeneous2D.cpp.

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

251  {
252  int cnt = 0, cnt1 = 0;
253  Array<OneD, NekDouble> tmparray;
254  int nlines = m_lines.num_elements();
255 
256  for(int n = 0; n < nlines; ++n)
257  {
258  m_lines[n]->FwdTrans(inarray+cnt, tmparray = outarray + cnt1,
259  coeffstate);
260  cnt += m_lines[n]->GetTotPoints();
261  cnt1 += m_lines[n]->GetNcoeffs();
262  }
263  if(!m_WaveSpace)
264  {
265  HomogeneousFwdTrans(outarray,outarray,coeffstate);
266  }
267  }
void HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
void Nektar::MultiRegions::ExpListHomogeneous2D::v_FwdTrans_IterPerExp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

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 269 of file ExpListHomogeneous2D.cpp.

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

270  {
271  int cnt = 0, cnt1 = 0;
272  Array<OneD, NekDouble> tmparray;
273  int nlines = m_lines.num_elements();
274 
275  for(int n = 0; n < nlines; ++n)
276  {
277  m_lines[n]->FwdTrans_IterPerExp(inarray+cnt, tmparray = outarray + cnt1);
278 
279  cnt += m_lines[n]->GetTotPoints();
280  cnt1 += m_lines[n]->GetNcoeffs();
281  }
282  if(!m_WaveSpace)
283  {
284  HomogeneousFwdTrans(outarray,outarray);
285  }
286  }
void HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
std::vector< LibUtilities::FieldDefinitionsSharedPtr > Nektar::MultiRegions::ExpListHomogeneous2D::v_GetFieldDefinitions ( void  )
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 582 of file ExpListHomogeneous2D.cpp.

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

583  {
584  std::vector<LibUtilities::FieldDefinitionsSharedPtr> returnval;
585  // Set up Homogeneous length details.
586  Array<OneD,LibUtilities::BasisSharedPtr> HomoBasis(2);
587  HomoBasis[0] = m_homogeneousBasis_y;
588  HomoBasis[1] = m_homogeneousBasis_z;
589 
590  std::vector<NekDouble> HomoLen(2);
591  HomoLen[0] = m_lhom_y;
592  HomoLen[1] = m_lhom_z;
593 
594  int nhom_modes_y = m_homogeneousBasis_y->GetNumModes();
595  int nhom_modes_z = m_homogeneousBasis_z->GetNumModes();
596 
597  std::vector<unsigned int> sIDs
599 
600  std::vector<unsigned int> yIDs;
601  std::vector<unsigned int> zIDs;
602 
603  for(int n = 0; n < nhom_modes_z; ++n)
604  {
605  for(int m = 0; m < nhom_modes_y; ++m)
606  {
607  zIDs.push_back(n);
608  yIDs.push_back(m);
609  }
610  }
611 
612  m_lines[0]->GeneralGetFieldDefinitions(returnval, 2, HomoBasis,
613  HomoLen, false,
614  sIDs, zIDs, yIDs);
615  return returnval;
616  }
NekDouble m_lhom_z
Width of homogeneous direction z.
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
LibUtilities::BasisSharedPtr m_homogeneousBasis_z
Base expansion in z direction.
NekDouble m_lhom_y
Width of homogeneous direction y.
static std::vector< unsigned int > NullUnsignedIntVector
Definition: FieldIO.h:51
LibUtilities::BasisSharedPtr m_homogeneousBasis_y
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
void Nektar::MultiRegions::ExpListHomogeneous2D::v_GetFieldDefinitions ( std::vector< LibUtilities::FieldDefinitionsSharedPtr > &  fielddef)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 618 of file ExpListHomogeneous2D.cpp.

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

619  {
620  // Set up Homogeneous length details.
621  Array<OneD,LibUtilities::BasisSharedPtr> HomoBasis(2);
622  HomoBasis[0] = m_homogeneousBasis_y;
623  HomoBasis[1] = m_homogeneousBasis_z;
624  std::vector<NekDouble> HomoLen(2);
625  HomoLen[0] = m_lhom_y;
626  HomoLen[1] = m_lhom_z;
627 
628  int nhom_modes_y = m_homogeneousBasis_y->GetNumModes();
629  int nhom_modes_z = m_homogeneousBasis_z->GetNumModes();
630 
631  std::vector<unsigned int> sIDs
633 
634  std::vector<unsigned int> yIDs;
635  std::vector<unsigned int> zIDs;
636 
637  for(int n = 0; n < nhom_modes_z; ++n)
638  {
639  for(int m = 0; m < nhom_modes_y; ++m)
640  {
641  zIDs.push_back(n);
642  yIDs.push_back(m);
643  }
644  }
645 
646  // enforce NumHomoDir == 1 by direct call
647  m_lines[0]->GeneralGetFieldDefinitions(fielddef, 2, HomoBasis,
648  HomoLen, false,
649  sIDs, zIDs, yIDs);
650  }
NekDouble m_lhom_z
Width of homogeneous direction z.
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
LibUtilities::BasisSharedPtr m_homogeneousBasis_z
Base expansion in z direction.
NekDouble m_lhom_y
Width of homogeneous direction y.
static std::vector< unsigned int > NullUnsignedIntVector
Definition: FieldIO.h:51
LibUtilities::BasisSharedPtr m_homogeneousBasis_y
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
virtual int Nektar::MultiRegions::ExpListHomogeneous2D::v_GetNumElmts ( void  )
inlineprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 168 of file ExpListHomogeneous2D.h.

169  {
170  return m_lines[0]->GetExpSize();
171  }
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
void Nektar::MultiRegions::ExpListHomogeneous2D::v_HomogeneousBwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
CoeffState  coeffstate = eLocal,
bool  Shuff = true,
bool  UnShuff = true 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 149 of file ExpListHomogeneous2D.cpp.

References Homogeneous2DTrans().

Referenced by HomogeneousBwdTrans().

154  {
155  // Backwards trans
156  Homogeneous2DTrans(inarray,outarray,false,coeffstate,Shuff,UnShuff);
157  }
void Homogeneous2DTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool IsForwards, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
void Nektar::MultiRegions::ExpListHomogeneous2D::v_HomogeneousFwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
CoeffState  coeffstate = eLocal,
bool  Shuff = true,
bool  UnShuff = true 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 139 of file ExpListHomogeneous2D.cpp.

References Homogeneous2DTrans().

Referenced by HomogeneousFwdTrans().

144  {
145  // Forwards trans
146  Homogeneous2DTrans(inarray,outarray,true,coeffstate,Shuff,UnShuff);
147  }
void Homogeneous2DTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool IsForwards, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
void Nektar::MultiRegions::ExpListHomogeneous2D::v_IProductWRTBase ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
CoeffState  coeffstate 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 327 of file ExpListHomogeneous2D.cpp.

References m_lines.

328  {
329  int cnt = 0, cnt1 = 0;
330  Array<OneD, NekDouble> tmparray;
331  int nlines = m_lines.num_elements();
332 
333  for(int n = 0; n < nlines; ++n)
334  {
335  m_lines[n]->IProductWRTBase(inarray+cnt, tmparray = outarray + cnt1,coeffstate);
336 
337  cnt += m_lines[n]->GetNcoeffs();
338  cnt1 += m_lines[n]->GetTotPoints();
339  }
340  }
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
void Nektar::MultiRegions::ExpListHomogeneous2D::v_IProductWRTBase_IterPerExp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

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 342 of file ExpListHomogeneous2D.cpp.

References m_lines.

343  {
344  int cnt = 0, cnt1 = 0;
345  Array<OneD, NekDouble> tmparray;
346  int nlines = m_lines.num_elements();
347 
348  for(int n = 0; n < nlines; ++n)
349  {
350  m_lines[n]->IProductWRTBase_IterPerExp(inarray+cnt, tmparray = outarray + cnt1);
351 
352  cnt += m_lines[n]->GetNcoeffs();
353  cnt1 += m_lines[n]->GetTotPoints();
354  }
355  }
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
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 
)
protectedvirtual

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 752 of file ExpListHomogeneous2D.cpp.

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

Referenced by PhysDeriv().

757  {
758  int nyzlines = m_lines.num_elements(); //number of Fourier points in the Fourier directions (nF_pts)
759  int npoints = inarray.num_elements(); //number of total points = n. of Fourier points * n. of points per line (nT_pts)
760  int n_points_line = npoints/nyzlines; //number of points per line
761 
762  Array<OneD, NekDouble> temparray(npoints);
763  Array<OneD, NekDouble> temparray1(npoints);
764  Array<OneD, NekDouble> temparray2(npoints);
765  Array<OneD, NekDouble> tmp1;
766  Array<OneD, NekDouble> tmp2;
767  Array<OneD, NekDouble> tmp3;
768 
769  for( int i=0 ; i<nyzlines ; i++ )
770  {
771  m_lines[i]->PhysDeriv( tmp1 = inarray + i*n_points_line ,tmp2 = out_d0 + i*n_points_line);
772  }
773 
775  {
776  if(m_WaveSpace)
777  {
778  temparray = inarray;
779  }
780  else
781  {
782  HomogeneousFwdTrans(inarray,temparray);
783  }
784  NekDouble sign = -1.0;
785  NekDouble beta;
786 
787  //along y
788  for(int i = 0; i < m_ny; i++)
789  {
790  beta = -sign*2*M_PI*(i/2)/m_lhom_y;
791 
792  for(int j = 0; j < m_nz; j++)
793  {
794  Vmath::Smul(n_points_line,beta,tmp1 = temparray + n_points_line*(i+j*m_ny),1, tmp2 = temparray1 + n_points_line*((i-int(sign))+j*m_ny),1);
795  }
796 
797  sign = -1.0*sign;
798  }
799 
800  //along z
801  sign = -1.0;
802  for(int i = 0; i < m_nz; i++)
803  {
804  beta = -sign*2*M_PI*(i/2)/m_lhom_z;
805  Vmath::Smul(m_ny*n_points_line,beta,tmp1 = temparray + i*m_ny*n_points_line,1,tmp2 = temparray2 + (i-int(sign))*m_ny*n_points_line,1);
806  sign = -1.0*sign;
807  }
808  if(m_WaveSpace)
809  {
810  out_d1 = temparray1;
811  out_d2 = temparray2;
812  }
813  else
814  {
815  HomogeneousBwdTrans(temparray1,out_d1);
816  HomogeneousBwdTrans(temparray2,out_d2);
817  }
818  }
819  else
820  {
821  if(m_WaveSpace)
822  {
823  ASSERTL0(false,"Semi-phyisical time-stepping not implemented yet for non-Fourier basis")
824  }
825  else
826  {
827  StdRegions::StdQuadExp StdQuad(m_homogeneousBasis_y->GetBasisKey(),m_homogeneousBasis_z->GetBasisKey());
828 
829  m_transposition->Transpose(inarray,temparray,false,LibUtilities::eXtoYZ);
830 
831  for(int i = 0; i < n_points_line; i++)
832  {
833  StdQuad.PhysDeriv(tmp1 = temparray + i*nyzlines, tmp2 = temparray1 + i*nyzlines, tmp3 = temparray2 + i*nyzlines);
834  }
835 
836  m_transposition->Transpose(temparray1,out_d1,false,LibUtilities::eYZtoX);
837  m_transposition->Transpose(temparray2,out_d2,false,LibUtilities::eYZtoX);
838  Vmath::Smul(npoints,2.0/m_lhom_y,out_d1,1,out_d1,1);
839  Vmath::Smul(npoints,2.0/m_lhom_z,out_d2,1,out_d2,1);
840  }
841  }
842  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
LibUtilities::TranspositionSharedPtr m_transposition
void HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:22
Fourier Expansion .
Definition: BasisType.h:52
NekDouble m_lhom_z
Width of homogeneous direction z.
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
LibUtilities::BasisSharedPtr m_homogeneousBasis_z
Base expansion in z direction.
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
NekDouble m_lhom_y
Width of homogeneous direction y.
double NekDouble
int m_ny
Number of modes = number of poitns in y direction.
LibUtilities::BasisSharedPtr m_homogeneousBasis_y
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
int m_nz
Number of modes = number of poitns in z direction.
void HomogeneousBwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
void Nektar::MultiRegions::ExpListHomogeneous2D::v_PhysDeriv ( Direction  edir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 844 of file ExpListHomogeneous2D.cpp.

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

848  {
849  int nyzlines = m_lines.num_elements(); //number of Fourier points in the Fourier directions (nF_pts)
850  int npoints = inarray.num_elements(); //number of total points = n. of Fourier points * n. of points per line (nT_pts)
851  int n_points_line = npoints/nyzlines; //number of points per line
852  //convert enum into int
853  int dir = (int)edir;
854 
855  Array<OneD, NekDouble> temparray(npoints);
856  Array<OneD, NekDouble> temparray1(npoints);
857  Array<OneD, NekDouble> temparray2(npoints);
858  Array<OneD, NekDouble> tmp1;
859  Array<OneD, NekDouble> tmp2;
860  Array<OneD, NekDouble> tmp3;
861 
862  if (dir < 1)
863  {
864  for( int i=0 ; i<nyzlines ; i++)
865  {
866  m_lines[i]->PhysDeriv( tmp1 = inarray + i*n_points_line ,tmp2 = out_d + i*n_points_line);
867  }
868  }
869  else
870  {
872  {
873  if(m_WaveSpace)
874  {
875  temparray = inarray;
876  }
877  else
878  {
879  HomogeneousFwdTrans(inarray,temparray);
880  }
881  NekDouble sign = -1.0;
882  NekDouble beta;
883 
884  if (dir == 1)
885  {
886  //along y
887  for(int i = 0; i < m_ny; i++)
888  {
889  beta = -sign*2*M_PI*(i/2)/m_lhom_y;
890 
891  for(int j = 0; j < m_nz; j++)
892  {
893  Vmath::Smul(n_points_line,beta,tmp1 = temparray + n_points_line*(i+j*m_ny),1, tmp2 = temparray1 + n_points_line*((i-int(sign))+j*m_ny),1);
894  }
895  sign = -1.0*sign;
896  }
897  if(m_WaveSpace)
898  {
899  out_d = temparray1;
900  }
901  else
902  {
903  HomogeneousBwdTrans(temparray1,out_d);
904  }
905  }
906  else
907  {
908  //along z
909  for(int i = 0; i < m_nz; i++)
910  {
911  beta = -sign*2*M_PI*(i/2)/m_lhom_z;
912  Vmath::Smul(m_ny*n_points_line,beta,tmp1 = temparray + i*m_ny*n_points_line,1,tmp2 = temparray2 + (i-int(sign))*m_ny*n_points_line,1);
913  sign = -1.0*sign;
914  }
915  if(m_WaveSpace)
916  {
917  out_d = temparray2;
918  }
919  else
920  {
921  HomogeneousBwdTrans(temparray2,out_d);
922  }
923  }
924  }
925  else
926  {
927  if(m_WaveSpace)
928  {
929  ASSERTL0(false,"Semi-phyisical time-stepping not implemented yet for non-Fourier basis")
930  }
931  else
932  {
933  StdRegions::StdQuadExp StdQuad(m_homogeneousBasis_y->GetBasisKey(),m_homogeneousBasis_z->GetBasisKey());
934 
935  m_transposition->Transpose(inarray,temparray,false,LibUtilities::eXtoYZ);
936 
937  for(int i = 0; i < n_points_line; i++)
938  {
939  StdQuad.PhysDeriv(tmp1 = temparray + i*nyzlines, tmp2 = temparray1 + i*nyzlines, tmp3 = temparray2 + i*nyzlines);
940  }
941 
942  if (dir == 1)
943  {
944  m_transposition->Transpose(temparray1,out_d,false,LibUtilities::eYZtoX);
945  Vmath::Smul(npoints,2.0/m_lhom_y,out_d,1,out_d,1);
946  }
947  else
948  {
949  m_transposition->Transpose(temparray2,out_d,false,LibUtilities::eYZtoX);
950  Vmath::Smul(npoints,2.0/m_lhom_z,out_d,1,out_d,1);
951  }
952  }
953  }
954  }
955  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
LibUtilities::TranspositionSharedPtr m_transposition
void HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:22
Fourier Expansion .
Definition: BasisType.h:52
NekDouble m_lhom_z
Width of homogeneous direction z.
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
LibUtilities::BasisSharedPtr m_homogeneousBasis_z
Base expansion in z direction.
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
NekDouble m_lhom_y
Width of homogeneous direction y.
double NekDouble
int m_ny
Number of modes = number of poitns in y direction.
LibUtilities::BasisSharedPtr m_homogeneousBasis_y
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
int m_nz
Number of modes = number of poitns in z direction.
void HomogeneousBwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
void Nektar::MultiRegions::ExpListHomogeneous2D::v_WriteVtkPieceData ( std::ostream &  outfile,
int  expansion,
std::string  var 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 729 of file ExpListHomogeneous2D.cpp.

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

731  {
732  int i;
733  int nq = (*m_exp)[expansion]->GetTotPoints();
734  int npoints_per_line = m_lines[0]->GetTotPoints();
735 
736  // printing the fields of that zone
737  outfile << " <DataArray type=\"Float64\" Name=\""
738  << var << "\">" << endl;
739  outfile << " ";
740  for (int n = 0; n < m_lines.num_elements(); ++n)
741  {
742  const Array<OneD, NekDouble> phys = m_phys + m_phys_offset[expansion] + n*npoints_per_line;
743  for(i = 0; i < nq; ++i)
744  {
745  outfile << (fabs(phys[i]) < NekConstants::kNekZeroTol ? 0 : phys[i]) << " ";
746  }
747  }
748  outfile << endl;
749  outfile << " </DataArray>" << endl;
750  }
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:956
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
static const NekDouble kNekZeroTol
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:991

Member Data Documentation

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

Definition at line 232 of file ExpListHomogeneous2D.h.

Referenced by ExpListHomogeneous2D().

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

Definition at line 138 of file ExpListHomogeneous2D.h.

Referenced by ExpListHomogeneous2D(), and Homogeneous2DTrans().

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

Definition at line 139 of file ExpListHomogeneous2D.h.

Referenced by ExpListHomogeneous2D(), and Homogeneous2DTrans().

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

Definition at line 158 of file ExpListHomogeneous2D.h.

Referenced by GetHomogeneous2DBlockMatrix().

LibUtilities::BasisSharedPtr Nektar::MultiRegions::ExpListHomogeneous2D::m_homogeneousBasis_y
protected
LibUtilities::BasisSharedPtr Nektar::MultiRegions::ExpListHomogeneous2D::m_homogeneousBasis_z
protected
NekDouble Nektar::MultiRegions::ExpListHomogeneous2D::m_lhom_y
protected
NekDouble Nektar::MultiRegions::ExpListHomogeneous2D::m_lhom_z
protected
Array<OneD, ExpListSharedPtr> Nektar::MultiRegions::ExpListHomogeneous2D::m_lines
protected

Vector of ExpList, will be filled with ExpList1D.

Definition at line 159 of file ExpListHomogeneous2D.h.

Referenced by Nektar::MultiRegions::ContField3DHomogeneous2D::ContField3DHomogeneous2D(), Nektar::MultiRegions::DisContField3DHomogeneous2D::DisContField3DHomogeneous2D(), Nektar::MultiRegions::DisContField3DHomogeneous2D::EvaluateBoundaryConditions(), Nektar::MultiRegions::ExpList1DHomogeneous2D::ExpList1DHomogeneous2D(), Nektar::MultiRegions::ExpList3DHomogeneous2D::ExpList3DHomogeneous2D(), ExpListHomogeneous2D(), GenHomogeneous2DBlockMatrix(), Nektar::MultiRegions::DisContField3DHomogeneous2D::GetBoundaryToElmtMap(), Nektar::MultiRegions::ExpList1DHomogeneous2D::GetCoords(), Homogeneous2DTrans(), Nektar::MultiRegions::ExpList1DHomogeneous2D::SetCoeffPhys(), Nektar::MultiRegions::ExpList3DHomogeneous2D::SetCoeffPhys(), Nektar::MultiRegions::DisContField3DHomogeneous2D::SetupBoundaryConditions(), v_AppendFieldData(), v_BwdTrans(), v_BwdTrans_IterPerExp(), Nektar::MultiRegions::ContField3DHomogeneous2D::v_ClearGlobalLinSysManager(), v_DealiasedProd(), v_ExtractDataToCoeffs(), v_FwdTrans(), v_FwdTrans_IterPerExp(), Nektar::MultiRegions::ExpList1DHomogeneous2D::v_GetCoords(), Nektar::MultiRegions::ExpList3DHomogeneous2D::v_GetCoords(), v_GetFieldDefinitions(), Nektar::MultiRegions::ContField3DHomogeneous2D::v_GlobalToLocal(), Nektar::MultiRegions::ContField3DHomogeneous2D::v_HelmSolve(), Nektar::MultiRegions::DisContField3DHomogeneous2D::v_HelmSolve(), Nektar::MultiRegions::ContField3DHomogeneous2D::v_ImposeDirichletConditions(), v_IProductWRTBase(), v_IProductWRTBase_IterPerExp(), Nektar::MultiRegions::ExpList3DHomogeneous2D::v_L2(), Nektar::MultiRegions::ContField3DHomogeneous2D::v_LocalToGlobal(), v_PhysDeriv(), and v_WriteVtkPieceData().

int Nektar::MultiRegions::ExpListHomogeneous2D::m_ny
protected
int Nektar::MultiRegions::ExpListHomogeneous2D::m_nz
protected
LibUtilities::BasisSharedPtr Nektar::MultiRegions::ExpListHomogeneous2D::m_paddingBasis_y
protected

Base expansion in y direction.

Definition at line 154 of file ExpListHomogeneous2D.h.

Referenced by SetPaddingBase().

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

Base expansion in z direction.

Definition at line 155 of file ExpListHomogeneous2D.h.

Referenced by SetPaddingBase().

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

Definition at line 233 of file ExpListHomogeneous2D.h.

Referenced by SetPaddingBase(), and v_DealiasedProd().

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

Definition at line 234 of file ExpListHomogeneous2D.h.

Referenced by SetPaddingBase(), and v_DealiasedProd().

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

Definition at line 140 of file ExpListHomogeneous2D.h.

Referenced by Homogeneous2DTrans().

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

Definition at line 141 of file ExpListHomogeneous2D.h.

Referenced by Homogeneous2DTrans().

LibUtilities::TranspositionSharedPtr Nektar::MultiRegions::ExpListHomogeneous2D::m_transposition
bool Nektar::MultiRegions::ExpListHomogeneous2D::m_useFFT
LibUtilities::CommSharedPtr Nektar::MultiRegions::ExpListHomogeneous2D::m_Ycomm

Definition at line 144 of file ExpListHomogeneous2D.h.

Referenced by ExpListHomogeneous2D().

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

Definition at line 145 of file ExpListHomogeneous2D.h.

Referenced by ExpListHomogeneous2D().

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

Definition at line 236 of file ExpListHomogeneous2D.h.

Referenced by SetPaddingBase().

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

Definition at line 235 of file ExpListHomogeneous2D.h.

Referenced by SetPaddingBase().