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::ExpListHomogeneous1D Class Reference

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

#include <ExpListHomogeneous1D.h>

Inheritance diagram for Nektar::MultiRegions::ExpListHomogeneous1D:
Inheritance graph
[legend]
Collaboration diagram for Nektar::MultiRegions::ExpListHomogeneous1D:
Collaboration graph
[legend]

Public Member Functions

 ExpListHomogeneous1D ()
 Default constructor.
 ExpListHomogeneous1D (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::BasisKey &HomoBasis, const NekDouble lz, const bool useFFT, const bool dealiasing)
 ExpListHomogeneous1D (const ExpListHomogeneous1D &In)
 Copy constructor.
virtual ~ExpListHomogeneous1D ()
 Destructor.
void Homogeneous1DTrans (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)
LibUtilities::BasisSharedPtr GetHomogeneousBasis (void)
void PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
 This function discretely evaluates the derivative of a function $f(\boldsymbol{x})$ on the domain consisting of all elements of the expansion.
void PhysDeriv (Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
ExpListSharedPtrGetPlane (int n)
- Public Member Functions inherited from Nektar::MultiRegions::ExpList
 ExpList ()
 The default constructor.
 ExpList (const LibUtilities::SessionReaderSharedPtr &pSession)
 The default constructor.
 ExpList (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 The default constructor.
 ExpList (const ExpList &in, const bool DeclareCoeffPhysArrays=true)
 The copy constructor.
virtual ~ExpList ()
 The default destructor.
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$.
int GetNcoeffs (const int eid) const
 Returns the total number of local degrees of freedom for element eid.
ExpansionType GetExpType (void)
 Returns the type of the expansion.
void SetExpType (ExpansionType Type)
 Returns the type of the expansion.
int EvalBasisNumModesMax (void) const
 Evaulates the maximum number of modes in the elemental basis order over all elements.
const Array< OneD, int > EvalBasisNumModesMaxPerExp (void) const
 Returns the vector of the number of modes in the elemental basis order over all elements.
int GetTotPoints (void) const
 Returns the total number of quadrature points m_npoints $=Q_{\mathrm{tot}}$.
int GetTotPoints (const int eid) const
 Returns the total number of quadrature points for eid's element $=Q_{\mathrm{tot}}$.
int GetNpoints (void) const
 Returns the total number of quadrature points m_npoints $=Q_{\mathrm{tot}}$.
int Get1DScaledTotPoints (const NekDouble scale) const
 Returns the total number of qudature points scaled by the factor scale on each 1D direction.
void SetWaveSpace (const bool wavespace)
 Sets the wave space to the one of the possible configuration true or false.
void SetModifiedBasis (const bool modbasis)
 Set Modified Basis for the stability analysis.
void SetPhys (int i, NekDouble val)
 Set the i th value of m_phys to value val.
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.
void SetPhys (const Array< OneD, const NekDouble > &inarray)
 Fills the array m_phys.
void SetPhysArray (Array< OneD, NekDouble > &inarray)
 Sets the array m_phys.
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.
bool GetPhysState (void) const
 This function indicates whether the array of physical values $\boldsymbol{u}_l$ (implemented as m_phys) is filled or not.
NekDouble PhysIntegral (void)
 This function integrates a function $f(\boldsymbol{x})$ over the domain consisting of all the elements of the expansion.
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.
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})$.
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.
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.
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.
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.
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.
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.
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.
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.
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$.
void GetBCValues (Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
void NormVectorIProductWRTBase (Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
void ApplyGeomInfo ()
 Apply geometry information to each expansion.
void WriteTecplotHeader (std::ofstream &outfile, std::string var="")
void WriteTecplotZone (std::ofstream &outfile, int expansion=-1)
void WriteTecplotField (std::ofstream &outfile, int expansion=-1)
void WriteTecplotConnectivity (std::ofstream &outfile, int expansion=-1)
void WriteVtkHeader (std::ofstream &outfile)
void WriteVtkFooter (std::ofstream &outfile)
void WriteVtkPieceHeader (std::ofstream &outfile, int expansion)
void WriteVtkPieceFooter (std::ofstream &outfile, int expansion)
void WriteVtkPieceData (std::ofstream &outfile, int expansion, std::string var="v")
int GetCoordim (int eid)
 This function returns the dimension of the coordinates of the element eid.
void SetCoeff (int i, NekDouble val)
 Set the i th coefficiient in m_coeffs to value val.
void SetCoeffs (int i, NekDouble val)
 Set the i th coefficiient in m_coeffs to value val.
void SetCoeffsArray (Array< OneD, NekDouble > &inarray)
 Set the m_coeffs array to inarray.
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.
void ImposeDirichletConditions (Array< OneD, NekDouble > &outarray)
 Impose Dirichlet Boundary Conditions onto Array.
void FillBndCondFromField (void)
 Fill Bnd Condition expansion from the values stored in expansion.
void LocalToGlobal (void)
 Put the coefficients into global ordering using m_coeffs.
void GlobalToLocal (void)
 Put the coefficients into local ordering and place in m_coeffs.
NekDouble GetCoeff (int i)
 Get the i th value (coefficient) of m_coeffs.
NekDouble GetCoeffs (int i)
 Get the i th value (coefficient) of m_coeffs.
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.
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.
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.
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.
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.
void SetHomo1DSpecVanVisc (Array< OneD, NekDouble > visc)
 This function sets the Spectral Vanishing Viscosity in homogeneous1D expansion.
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.
LibUtilities::TranspositionSharedPtr GetTransposition (void)
 This function returns the transposition class associaed with the homogeneous expansion.
NekDouble GetHomoLen (void)
 This function returns the Width of homogeneous direction associaed with the homogeneous expansion.
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.
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.
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.
int GetExpSize (void)
 This function returns the number of elements in the expansion.
int GetNumElmts (void)
 This function returns the number of elements in the expansion which may be different for a homogeoenous extended expansionp.
const boost::shared_ptr
< LocalRegions::ExpansionVector
GetExp () const
 This function returns the vector of elements in the expansion.
LocalRegions::ExpansionSharedPtrGetExp (int n) const
 This function returns (a shared pointer to) the local elemental expansion of the $n^{\mathrm{th}}$ element.
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.
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.
int GetPhys_Offset (int n) const
 Get the start offset position for a global list of m_phys correspoinding to element n.
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.
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.
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.
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
const Array< OneD, const
boost::shared_ptr< ExpList > > & 
GetBndCondExpansions ()
boost::shared_ptr< ExpList > & UpdateBndCondExpansion (int i)
void Upwind (const Array< OneD, const Array< OneD, NekDouble > > &Vec, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
void Upwind (const Array< OneD, const NekDouble > &Vn, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
boost::shared_ptr< ExpList > & GetTrace ()
boost::shared_ptr
< AssemblyMapDG > & 
GetTraceMap (void)
const Array< OneD, const int > & GetTraceBndMap (void)
void GetNormals (Array< OneD, Array< OneD, NekDouble > > &normals)
void AddTraceIntegral (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
void AddTraceIntegral (const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
void AddFwdBwdTraceIntegral (const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
void GetFwdBwdTracePhys (Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
void GetFwdBwdTracePhys (const Array< OneD, const NekDouble > &field, Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
void ExtractTracePhys (Array< OneD, NekDouble > &outarray)
void ExtractTracePhys (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
const Array< OneD, const
SpatialDomains::BoundaryConditionShPtr > & 
GetBndConditions ()
Array< OneD,
SpatialDomains::BoundaryConditionShPtr > & 
UpdateBndConditions ()
void EvaluateBoundaryConditions (const NekDouble time=0.0, const std::string varName="", const NekDouble=NekConstants::kNekUnsetDouble, const NekDouble=NekConstants::kNekUnsetDouble)
void GeneralMatrixOp (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
 This function calculates the result of the multiplication of a matrix of type specified by mkey with a vector given by inarray.
void GeneralMatrixOp_IterPerExp (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void SetUpPhysNormals ()
void GetBoundaryToElmtMap (Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
void GeneralGetFieldDefinitions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef, int NumHomoDir=0, Array< OneD, LibUtilities::BasisSharedPtr > &HomoBasis=LibUtilities::NullBasisSharedPtr1DArray, std::vector< NekDouble > &HomoLen=LibUtilities::NullNekDoubleVector, std::vector< unsigned int > &HomoZIDs=LibUtilities::NullUnsignedIntVector, std::vector< unsigned int > &HomoYIDs=LibUtilities::NullUnsignedIntVector)
const
NekOptimize::GlobalOptParamSharedPtr
GetGlobalOptParam (void)
map< int, RobinBCInfoSharedPtrGetRobinBCInfo ()
void GetPeriodicEntities (PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces=NullPeriodicMap)
std::vector
< LibUtilities::FieldDefinitionsSharedPtr
GetFieldDefinitions ()
void GetFieldDefinitions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef)
void AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
 Append the element data listed in elements fielddef->m_ElementIDs onto fielddata.
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.
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.
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.
void ExtractDataToCoeffs (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
 Extract the data in fielddata into the coeffs.
boost::shared_ptr< ExpListGetSharedThisPtr ()
 Returns a shared pointer to the current object.
boost::shared_ptr
< LibUtilities::SessionReader
GetSession ()
 Returns the session object.
boost::shared_ptr
< LibUtilities::Comm
GetComm ()
 Returns the comm object.
SpatialDomains::MeshGraphSharedPtr GetGraph ()

Public Attributes

LibUtilities::TranspositionSharedPtr m_transposition
- Public Attributes inherited from Nektar::MultiRegions::ExpList
ExpansionType m_expType

Protected Member Functions

DNekBlkMatSharedPtr GenHomogeneous1DBlockMatrix (Homogeneous1DMatType mattype, CoeffState coeffstate=eLocal) const
DNekBlkMatSharedPtr GetHomogeneous1DBlockMatrix (Homogeneous1DMatType mattype, CoeffState coeffstate=eLocal) const
NekDouble GetSpecVanVisc (const int k)
virtual void v_SetHomo1DSpecVanVisc (Array< OneD, NekDouble > visc)
virtual int v_GetNumElmts (void)
virtual
LibUtilities::BasisSharedPtr 
v_GetHomogeneousBasis (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.
virtual void v_ExtractCoeffsToCoeffs (const boost::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
virtual void v_WriteVtkPieceData (std::ofstream &outfile, int expansion, std::string var)
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_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)
virtual
LibUtilities::TranspositionSharedPtr 
v_GetTransposition (void)
virtual Array< OneD, const
unsigned int > 
v_GetZIDs (void)
virtual ExpListSharedPtrv_GetPlane (int n)
virtual NekDouble v_GetHomoLen (void)
- 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.
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.
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.
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.
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 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_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_SetUpPhysNormals ()
virtual void v_GetBoundaryToElmtMap (Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
virtual void v_ReadGlobalOptimizationParameters ()
virtual void v_WriteTecplotHeader (std::ofstream &outfile, std::string var="")
virtual void v_WriteTecplotZone (std::ofstream &outfile, int expansion)
virtual void v_WriteTecplotField (std::ofstream &outfile, int expansion)
virtual void v_WriteTecplotConnectivity (std::ofstream &outfile, int expansion)
virtual void v_WriteVtkPieceHeader (std::ofstream &outfile, int expansion)
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 Array< OneD, const
unsigned int > 
v_GetYIDs (void)
void ExtractFileBCs (const std::string &fileName, const std::string &varName, const boost::shared_ptr< ExpList > locExpList)

Protected Attributes

bool m_useFFT
 FFT variables.
LibUtilities::NektarFFTSharedPtr m_FFT
LibUtilities::NektarFFTSharedPtr m_FFT_deal
Array< OneD, NekDoublem_tmpIN
Array< OneD, NekDoublem_tmpOUT
LibUtilities::BasisSharedPtr m_homogeneousBasis
 Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m_coeff and m_phys.
NekDouble m_lhom
 Width of homogeneous direction.
Homo1DBlockMatrixMapShPtr m_homogeneous1DBlockMat
Array< OneD, ExpListSharedPtrm_planes
- Protected Attributes inherited from Nektar::MultiRegions::ExpList
LibUtilities::CommSharedPtr m_comm
 Communicator.
LibUtilities::SessionReaderSharedPtr m_session
 Session.
SpatialDomains::MeshGraphSharedPtr m_graph
 Mesh associated with this expansion list.
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$.
int m_npoints
Array< OneD, NekDoublem_coeffs
 Concatenation of all local expansion coefficients.
Array< OneD, NekDoublem_phys
 The global expansion evaluated at the quadrature points.
bool m_physState
 The state of the array m_phys.
boost::shared_ptr
< LocalRegions::ExpansionVector
m_exp
 The list of local expansions.
Array< OneD, int > m_coeff_offset
 Offset of elemental data into the array m_coeffs.
Array< OneD, int > m_phys_offset
 Offset of elemental data into the array m_phys.
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]];.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
BlockMatrixMapShPtr m_blockMat
bool m_WaveSpace

Private Attributes

bool m_dealiasing
int m_padsize
Array< OneD, NekDoublem_specVanVisc
 Spectral vanishing Viscosity coefficient for stabilisation.

Additional Inherited Members

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

Detailed Description

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

Definition at line 77 of file ExpListHomogeneous1D.h.

Constructor & Destructor Documentation

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

Definition at line 56 of file ExpListHomogeneous1D.cpp.

References ASSERTL0, ASSERTL2, Nektar::LibUtilities::BasisManager(), Nektar::LibUtilities::NekFactory< tKey, tBase, >::CreateInstance(), Nektar::LibUtilities::GetNektarFFTFactory(), Nektar::MultiRegions::ExpList::m_comm, m_dealiasing, m_FFT, m_FFT_deal, m_homogeneousBasis, m_padsize, m_planes, m_transposition, m_useFFT, and Nektar::LibUtilities::NullBasisKey().

:
ExpList(pSession),
m_useFFT(useFFT),
m_lhom(lhom),
m_dealiasing(dealiasing)
{
ASSERTL2(HomoBasis != LibUtilities::NullBasisKey,"Homogeneous Basis is a null basis");
m_planes = Array<OneD,ExpListSharedPtr>(m_homogeneousBasis->GetNumPoints()/m_comm->GetColumnComm()->GetSize());
{
}
{
{
NekDouble size = 1.5*m_homogeneousBasis->GetNumPoints();
m_padsize = int(size);
.CreateInstance("NekFFTW", m_padsize);
}
else
{
ASSERTL0(false, "Dealiasing available just in combination "
"with FFTW");
}
}
}
Nektar::MultiRegions::ExpListHomogeneous1D::ExpListHomogeneous1D ( const ExpListHomogeneous1D In)

Copy constructor.

Parameters
InExpListHomogeneous1D object to copy.

Definition at line 98 of file ExpListHomogeneous1D.cpp.

References m_planes.

:
ExpList(In,false),
m_transposition(In.m_transposition),
m_useFFT(In.m_useFFT),
m_FFT(In.m_FFT),
m_tmpIN(In.m_tmpIN),
m_tmpOUT(In.m_tmpOUT),
m_homogeneousBasis(In.m_homogeneousBasis),
m_lhom(In.m_lhom),
m_homogeneous1DBlockMat(In.m_homogeneous1DBlockMat),
m_dealiasing(In.m_dealiasing),
m_padsize(In.m_padsize)
{
m_planes = Array<OneD, ExpListSharedPtr>(In.m_planes.num_elements());
}
Nektar::MultiRegions::ExpListHomogeneous1D::~ExpListHomogeneous1D ( )
virtual

Destructor.

Destructor

Definition at line 117 of file ExpListHomogeneous1D.cpp.

{
}

Member Function Documentation

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

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 297 of file ExpListHomogeneous1D.h.

References v_DealiasedProd().

{
v_DealiasedProd(inarray1,inarray2,outarray,coeffstate);
}
DNekBlkMatSharedPtr Nektar::MultiRegions::ExpListHomogeneous1D::GenHomogeneous1DBlockMatrix ( Homogeneous1DMatType  mattype,
CoeffState  coeffstate = eLocal 
) const
protected

Definition at line 479 of file ExpListHomogeneous1D.cpp.

References Nektar::MultiRegions::eBackwardsCoeffSpace1D, Nektar::StdRegions::eBwdTrans, Nektar::eDIAGONAL, Nektar::MultiRegions::eForwardsCoeffSpace1D, Nektar::MultiRegions::eForwardsPhysSpace1D, Nektar::LibUtilities::eFourierHalfModeIm, Nektar::LibUtilities::eFourierHalfModeRe, Nektar::StdRegions::eFwdTrans, Nektar::MultiRegions::ExpList::m_comm, m_homogeneousBasis, and m_planes.

Referenced by GetHomogeneous1DBlockMatrix().

{
int n_exp = 0;
int num_trans_per_proc = 0;
if((mattype == eForwardsCoeffSpace1D)
||(mattype == eBackwardsCoeffSpace1D)) // will operate on m_coeffs
{
n_exp = m_planes[0]->GetNcoeffs();
}
else
{
n_exp = m_planes[0]->GetTotPoints(); // will operatore on m_phys
}
num_trans_per_proc = n_exp/m_comm->GetColumnComm()->GetSize() + (n_exp%m_comm->GetColumnComm()->GetSize() > 0);
Array<OneD,unsigned int> nrows(num_trans_per_proc);
Array<OneD,unsigned int> ncols(num_trans_per_proc);
if((mattype == eForwardsCoeffSpace1D)||(mattype == eForwardsPhysSpace1D))
{
nrows = Array<OneD, unsigned int>(num_trans_per_proc,m_homogeneousBasis->GetNumModes());
ncols = Array<OneD, unsigned int>(num_trans_per_proc,m_homogeneousBasis->GetNumPoints());
}
else
{
nrows = Array<OneD, unsigned int>(num_trans_per_proc,m_homogeneousBasis->GetNumPoints());
ncols = Array<OneD, unsigned int>(num_trans_per_proc,m_homogeneousBasis->GetNumModes());
}
MatrixStorage blkmatStorage = eDIAGONAL;
::AllocateSharedPtr(nrows,ncols,blkmatStorage);
//Half Mode
{
StdRegions::StdPointExp StdPoint(m_homogeneousBasis->GetBasisKey());
if((mattype == eForwardsCoeffSpace1D)||(mattype == eForwardsPhysSpace1D))
{
StdPoint.DetShapeType(),
StdPoint);
loc_mat = StdPoint.GetStdMatrix(matkey);
}
else
{
StdPoint.DetShapeType(),
StdPoint);
loc_mat = StdPoint.GetStdMatrix(matkey);
}
}
//other cases
else
{
StdRegions::StdSegExp StdSeg(m_homogeneousBasis->GetBasisKey());
if((mattype == eForwardsCoeffSpace1D)||(mattype == eForwardsPhysSpace1D))
{
StdSeg.DetShapeType(),
StdSeg);
loc_mat = StdSeg.GetStdMatrix(matkey);
}
else
{
StdSeg.DetShapeType(),
StdSeg);
loc_mat = StdSeg.GetStdMatrix(matkey);
}
}
// set up array of block matrices.
for(int i = 0; i < num_trans_per_proc; ++i)
{
BlkMatrix->SetBlock(i,i,loc_mat);
}
return BlkMatrix;
}
DNekBlkMatSharedPtr Nektar::MultiRegions::ExpListHomogeneous1D::GetHomogeneous1DBlockMatrix ( Homogeneous1DMatType  mattype,
CoeffState  coeffstate = eLocal 
) const
protected

Definition at line 463 of file ExpListHomogeneous1D.cpp.

References GenHomogeneous1DBlockMatrix(), Nektar::iterator, and m_homogeneous1DBlockMat.

Referenced by Homogeneous1DTrans().

{
if(matrixIter == m_homogeneous1DBlockMat->end())
{
return ((*m_homogeneous1DBlockMat)[mattype] =
GenHomogeneous1DBlockMatrix(mattype,coeffstate));
}
else
{
return matrixIter->second;
}
}
LibUtilities::BasisSharedPtr Nektar::MultiRegions::ExpListHomogeneous1D::GetHomogeneousBasis ( void  )
inline

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 117 of file ExpListHomogeneous1D.h.

References m_homogeneousBasis.

Referenced by v_GetHomogeneousBasis().

{
}
ExpListSharedPtr& Nektar::MultiRegions::ExpListHomogeneous1D::GetPlane ( int  n)
inline

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 131 of file ExpListHomogeneous1D.h.

References m_planes.

Referenced by v_GetPlane().

{
return m_planes[n];
}
NekDouble Nektar::MultiRegions::ExpListHomogeneous1D::GetSpecVanVisc ( const int  k)
inlineprotected

Definition at line 162 of file ExpListHomogeneous1D.h.

References m_specVanVisc.

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

{
NekDouble returnval = 0.0;
if(m_specVanVisc.num_elements())
{
returnval = m_specVanVisc[k];
}
return returnval;
}
void Nektar::MultiRegions::ExpListHomogeneous1D::Homogeneous1DTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
bool  IsForwards,
CoeffState  coeffstate = eLocal,
bool  Shuff = true,
bool  UnShuff = true 
)

Homogeneous transform Bwd/Fwd (MVM and FFT)

Definition at line 341 of file ExpListHomogeneous1D.cpp.

References Nektar::MultiRegions::eBackwardsCoeffSpace1D, Nektar::MultiRegions::eBackwardsPhysSpace1D, Nektar::MultiRegions::eForwardsCoeffSpace1D, Nektar::MultiRegions::eForwardsPhysSpace1D, Nektar::eWrapper, Nektar::LibUtilities::eXYtoZ, Nektar::LibUtilities::eZtoXY, GetHomogeneous1DBlockMatrix(), Nektar::MultiRegions::ExpList::m_comm, m_FFT, m_homogeneousBasis, Nektar::MultiRegions::ExpList::m_npoints, m_planes, m_tmpIN, m_tmpOUT, m_transposition, m_useFFT, and Vmath::Vcopy().

Referenced by v_HomogeneousBwdTrans(), and v_HomogeneousFwdTrans().

{
int num_dofs;
if(IsForwards)
{
num_dofs = inarray.num_elements();
}
else
{
num_dofs = outarray.num_elements();
}
{
int num_points_per_plane = num_dofs/m_planes.num_elements();
int num_dfts_per_proc = num_points_per_plane/m_comm->GetColumnComm()->GetSize() + (num_points_per_plane%m_comm->GetColumnComm()->GetSize() > 0);
Array<OneD, NekDouble> fft_in (num_dfts_per_proc*m_homogeneousBasis->GetNumPoints(),0.0);
Array<OneD, NekDouble> fft_out(num_dfts_per_proc*m_homogeneousBasis->GetNumPoints(),0.0);
if(Shuff)
{
m_transposition->Transpose(inarray,fft_in,false,LibUtilities::eXYtoZ);
}
else
{
Vmath::Vcopy(num_dfts_per_proc*m_homogeneousBasis->GetNumPoints(),
inarray,1,fft_in,1);
}
if(IsForwards)
{
for(int i = 0 ; i < num_dfts_per_proc ; i++)
{
m_FFT->FFTFwdTrans(m_tmpIN = fft_in + i*m_homogeneousBasis->GetNumPoints(), m_tmpOUT = fft_out + i*m_homogeneousBasis->GetNumPoints());
}
}
else
{
for(int i = 0 ; i < num_dfts_per_proc ; i++)
{
m_FFT->FFTBwdTrans(m_tmpIN = fft_in + i*m_homogeneousBasis->GetNumPoints(), m_tmpOUT = fft_out + i*m_homogeneousBasis->GetNumPoints());
}
}
if(UnShuff)
{
m_transposition->Transpose(fft_out,outarray,false,LibUtilities::eZtoXY);
}
else
{
Vmath::Vcopy(num_dfts_per_proc*m_homogeneousBasis->GetNumPoints(),
fft_out,1,outarray,1);
}
}
else
{
if(num_dofs == m_npoints) //transform phys space
{
if(IsForwards)
{
}
else
{
}
}
else
{
if(IsForwards)
{
}
else
{
}
}
int nrows = blkmat->GetRows();
int ncols = blkmat->GetColumns();
Array<OneD, NekDouble> sortedinarray(ncols,0.0);
Array<OneD, NekDouble> sortedoutarray(nrows,0.0);
if(Shuff)
{
m_transposition->Transpose(inarray,sortedinarray,!IsForwards,LibUtilities::eXYtoZ);
}
else
{
Vmath::Vcopy(ncols,inarray,1,sortedinarray,1);
}
// Create NekVectors from the given data arrays
NekVector<NekDouble> in (ncols,sortedinarray,eWrapper);
NekVector<NekDouble> out(nrows,sortedoutarray,eWrapper);
// Perform matrix-vector multiply.
out = (*blkmat)*in;
if(UnShuff)
{
m_transposition->Transpose(sortedoutarray,outarray,IsForwards,LibUtilities::eZtoXY);
}
else
{
Vmath::Vcopy(nrows,sortedinarray,1,outarray,1);
}
}
}
void Nektar::MultiRegions::ExpListHomogeneous1D::HomogeneousBwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
CoeffState  coeffstate = eLocal,
bool  Shuff = true,
bool  UnShuff = true 
)
inline

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 288 of file ExpListHomogeneous1D.h.

References v_HomogeneousBwdTrans().

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

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

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 279 of file ExpListHomogeneous1D.h.

References v_HomogeneousFwdTrans().

Referenced by v_DealiasedProd(), v_FwdTrans(), v_FwdTrans_IterPerExp(), Nektar::MultiRegions::ContField3DHomogeneous1D::v_HelmSolve(), Nektar::MultiRegions::DisContField3DHomogeneous1D::v_HelmSolve(), and v_PhysDeriv().

{
v_HomogeneousFwdTrans(inarray,outarray,coeffstate,Shuff,UnShuff);
}
void Nektar::MultiRegions::ExpListHomogeneous1D::PhysDeriv ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d0,
Array< OneD, NekDouble > &  out_d1,
Array< OneD, NekDouble > &  out_d2 
)

This function discretely evaluates the derivative of a function $f(\boldsymbol{x})$ on the domain consisting of all elements of the expansion.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1047 of file ExpListHomogeneous1D.cpp.

References v_PhysDeriv().

{
v_PhysDeriv(inarray,out_d0,out_d1,out_d2);
}
void Nektar::MultiRegions::ExpListHomogeneous1D::PhysDeriv ( Direction  edir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d 
)

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1056 of file ExpListHomogeneous1D.cpp.

References v_PhysDeriv().

{
v_PhysDeriv(edir,inarray,out_d);
}
void Nektar::MultiRegions::ExpListHomogeneous1D::v_AppendFieldData ( LibUtilities::FieldDefinitionsSharedPtr fielddef,
std::vector< NekDouble > &  fielddata 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 645 of file ExpListHomogeneous1D.cpp.

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

{
v_AppendFieldData(fielddef,fielddata,m_coeffs);
}
void Nektar::MultiRegions::ExpListHomogeneous1D::v_AppendFieldData ( LibUtilities::FieldDefinitionsSharedPtr fielddef,
std::vector< NekDouble > &  fielddata,
Array< OneD, NekDouble > &  coeffs 
)
protectedvirtual

This routine appends the data from the expansion list into the output format where each element is given by looping over its Fourier modes where as data in the expandion is stored with all consecutive elements and then the Fourier modes

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 620 of file ExpListHomogeneous1D.cpp.

References Nektar::MultiRegions::ExpList::m_coeff_offset, and m_planes.

{
int i,n;
int ncoeffs_per_plane = m_planes[0]->GetNcoeffs();
// Determine mapping from element ids to location in
// expansion list
map<int, int> ElmtID_to_ExpID;
for(i = 0; i < m_planes[0]->GetExpSize(); ++i)
{
ElmtID_to_ExpID[(*m_exp)[i]->GetGeom()->GetGlobalID()] = i;
}
for(i = 0; i < fielddef->m_elementIDs.size(); ++i)
{
int eid = ElmtID_to_ExpID[fielddef->m_elementIDs[i]];
int datalen = (*m_exp)[eid]->GetNcoeffs();
for(n = 0; n < m_planes.num_elements(); ++n)
{
fielddata.insert(fielddata.end(),&coeffs[m_coeff_offset[eid]+n*ncoeffs_per_plane],&coeffs[m_coeff_offset[eid]+n*ncoeffs_per_plane]+datalen);
}
}
}
void Nektar::MultiRegions::ExpListHomogeneous1D::v_BwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
CoeffState  coeffstate 
)
protectedvirtual

Backward transform

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 265 of file ExpListHomogeneous1D.cpp.

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

{
int cnt = 0, cnt1 = 0;
Array<OneD, NekDouble> tmparray;
for(int n = 0; n < m_planes.num_elements(); ++n)
{
m_planes[n]->BwdTrans(inarray+cnt, tmparray = outarray + cnt1,
coeffstate);
cnt += m_planes[n]->GetNcoeffs();
cnt1 += m_planes[n]->GetTotPoints();
}
{
HomogeneousBwdTrans(outarray,outarray);
}
}
void Nektar::MultiRegions::ExpListHomogeneous1D::v_BwdTrans_IterPerExp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Backward transform element by element

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 286 of file ExpListHomogeneous1D.cpp.

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

{
int cnt = 0, cnt1 = 0;
Array<OneD, NekDouble> tmparray;
for(int n = 0; n < m_planes.num_elements(); ++n)
{
m_planes[n]->BwdTrans_IterPerExp(inarray+cnt, tmparray = outarray + cnt1);
cnt += m_planes[n]->GetNcoeffs();
cnt1 += m_planes[n]->GetTotPoints();
}
{
HomogeneousBwdTrans(outarray,outarray);
}
}
void Nektar::MultiRegions::ExpListHomogeneous1D::v_DealiasedProd ( const Array< OneD, NekDouble > &  inarray1,
const Array< OneD, NekDouble > &  inarray2,
Array< OneD, NekDouble > &  outarray,
CoeffState  coeffstate = eLocal 
)
protectedvirtual

Dealiasing routine

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 144 of file ExpListHomogeneous1D.cpp.

References Nektar::LibUtilities::eXYtoZ, Nektar::LibUtilities::eZtoXY, HomogeneousBwdTrans(), HomogeneousFwdTrans(), Nektar::MultiRegions::ExpList::m_comm, m_FFT_deal, m_homogeneousBasis, m_padsize, m_planes, m_transposition, Vmath::Vcopy(), and Vmath::Vmul().

Referenced by DealiasedProd().

{
// inarray1 = first term of the product in full physical space
// inarray2 = second term of the product in full physical space
// dealiased product stored in outarray
int num_dofs = inarray1.num_elements();
int N = m_homogeneousBasis->GetNumPoints();
Array<OneD, NekDouble> V1(num_dofs);
Array<OneD, NekDouble> V2(num_dofs);
Array<OneD, NekDouble> V1V2(num_dofs);
HomogeneousFwdTrans(inarray1,V1,coeffstate);
HomogeneousFwdTrans(inarray2,V2,coeffstate);
int num_points_per_plane = num_dofs/m_planes.num_elements();
int num_proc = m_comm->GetColumnComm()->GetSize();
int num_dfts_per_proc = num_points_per_plane / num_proc
+ (num_points_per_plane % num_proc > 0);
Array<OneD, NekDouble> ShufV1(num_dfts_per_proc*N,0.0);
Array<OneD, NekDouble> ShufV2(num_dfts_per_proc*N,0.0);
Array<OneD, NekDouble> ShufV1V2(num_dfts_per_proc*N,0.0);
Array<OneD, NekDouble> ShufV1_PAD_coef(m_padsize,0.0);
Array<OneD, NekDouble> ShufV2_PAD_coef(m_padsize,0.0);
Array<OneD, NekDouble> ShufV1_PAD_phys(m_padsize,0.0);
Array<OneD, NekDouble> ShufV2_PAD_phys(m_padsize,0.0);
Array<OneD, NekDouble> ShufV1V2_PAD_coef(m_padsize,0.0);
Array<OneD, NekDouble> ShufV1V2_PAD_phys(m_padsize,0.0);
m_transposition->Transpose(V1, ShufV1, false, LibUtilities::eXYtoZ);
m_transposition->Transpose(V2, ShufV2, false, LibUtilities::eXYtoZ);
// Looping on the pencils
for(int i = 0 ; i < num_dfts_per_proc ; i++)
{
// Copying the i-th pencil pf lenght N into a bigger
// pencil of lenght 2N We are in Fourier space
Vmath::Vcopy(N, &(ShufV1[i*N]), 1, &(ShufV1_PAD_coef[0]), 1);
Vmath::Vcopy(N, &(ShufV2[i*N]), 1, &(ShufV2_PAD_coef[0]), 1);
// Moving to physical space using the padded system
m_FFT_deal->FFTBwdTrans(ShufV1_PAD_coef, ShufV1_PAD_phys);
m_FFT_deal->FFTBwdTrans(ShufV2_PAD_coef, ShufV2_PAD_phys);
// Perfroming the vectors multiplication in physical space on
// the padded system
Vmath::Vmul(m_padsize, ShufV1_PAD_phys, 1,
ShufV2_PAD_phys, 1,
ShufV1V2_PAD_phys, 1);
// Moving back the result (V1*V2)_phys in Fourier space, padded
// system
m_FFT_deal->FFTFwdTrans(ShufV1V2_PAD_phys, ShufV1V2_PAD_coef);
// Copying the first half of the padded pencil in the full
// vector (Fourier space)
Vmath::Vcopy(N, &(ShufV1V2_PAD_coef[0]), 1,
&(ShufV1V2[i*N]), 1);
}
m_transposition->Transpose(ShufV1V2, V1V2, false,
// Moving the results in physical space for the output
HomogeneousBwdTrans(V1V2, outarray, coeffstate);
}
void Nektar::MultiRegions::ExpListHomogeneous1D::v_ExtractCoeffsToCoeffs ( const boost::shared_ptr< ExpList > &  fromExpList,
const Array< OneD, const NekDouble > &  fromCoeffs,
Array< OneD, NekDouble > &  toCoeffs 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 776 of file ExpListHomogeneous1D.cpp.

References Nektar::MultiRegions::ExpList::m_ncoeffs, and m_planes.

{
int i;
int fromNcoeffs_per_plane = fromExpList->GetPlane(0)->GetNcoeffs();
Array<OneD, NekDouble> tocoeffs_tmp, fromcoeffs_tmp;
for(i = 0; i < m_planes.num_elements(); ++i)
{
m_planes[i]->ExtractCoeffsToCoeffs(fromExpList->GetPlane(i),fromcoeffs_tmp = fromCoeffs + fromNcoeffs_per_plane*i, tocoeffs_tmp = toCoeffs + m_ncoeffs*i);
}
}
void Nektar::MultiRegions::ExpListHomogeneous1D::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 651 of file ExpListHomogeneous1D.cpp.

References ASSERTL1, Nektar::LibUtilities::GetNumberOfCoefficients(), Nektar::iterator, Nektar::MultiRegions::ExpList::m_coeff_offset, Nektar::MultiRegions::ExpList::m_exp, m_homogeneousBasis, m_planes, m_transposition, and Vmath::Vcopy().

{
int i,n;
int offset = 0;
int nzmodes;
int datalen = fielddata.size()/fielddef->m_fields.size();
std::vector<unsigned int> fieldDefHomoZids;
// Find data location according to field definition
for(i = 0; i < fielddef->m_fields.size(); ++i)
{
if(fielddef->m_fields[i] == field)
{
break;
}
offset += datalen;
}
if(i == fielddef->m_fields.size())
{
cout << "Field "<< field<< "not found in data file. " << endl;
}
else
{
int modes_offset = 0;
int planes_offset = 0;
Array<OneD, NekDouble> coeff_tmp;
// Build map of plane IDs lying on this processor.
std::map<int,int> homoZids;
for (i = 0; i < m_planes.num_elements(); ++i)
{
homoZids[m_transposition->GetPlaneID(i)] = i;
}
if(fielddef->m_numHomogeneousDir)
{
for(i = 0; i < fielddef->m_basis.size(); ++i)
{
if(fielddef->m_basis[i] == m_homogeneousBasis->GetBasisType())
{
nzmodes = fielddef->m_homogeneousZIDs.size();
break;
}
}
ASSERTL1(i != fielddef->m_basis.size(),"Failed to determine number of Homogeneous modes");
fieldDefHomoZids = fielddef->m_homogeneousZIDs;
}
else // input file is 2D and so set nzmodes to 1
{
nzmodes = 1;
fieldDefHomoZids.push_back(0);
}
// Determine mapping from element ids to location in expansion list.
map<int, int> ElmtID_to_ExpID;
for(i = 0; i < m_planes[0]->GetExpSize(); ++i)
{
ElmtID_to_ExpID[(*m_exp)[i]->GetGeom()->GetGlobalID()] = i;
}
// calculate number of modes in the current partition
int ncoeffs_per_plane = m_planes[0]->GetNcoeffs();
for(i = 0; i < fielddef->m_elementIDs.size(); ++i)
{
if(fielddef->m_uniOrder == true) // reset modes_offset to zero
{
modes_offset = 0;
}
int datalen = LibUtilities::GetNumberOfCoefficients(fielddef->m_shapeType,
fielddef->m_numModes,
modes_offset);
it = ElmtID_to_ExpID.find(fielddef->m_elementIDs[i]);
// ensure element is on this partition for parallel case.
if(it == ElmtID_to_ExpID.end())
{
// increase offset for correct FieldData access
offset += datalen*nzmodes;
continue;
}
int eid = it->second;
for(n = 0; n < nzmodes; ++n, offset += datalen)
{
it = homoZids.find(fieldDefHomoZids[n]);
// Check to make sure this mode number lies in this field.
if (it == homoZids.end())
{
continue;
}
planes_offset = it->second;
if(datalen == (*m_exp)[eid]->GetNcoeffs())
{
Vmath::Vcopy(datalen,&fielddata[offset],1,&coeffs[m_coeff_offset[eid]+planes_offset*ncoeffs_per_plane],1);
}
else // unpack data to new order
{
(*m_exp)[eid]->ExtractDataToCoeffs(&fielddata[offset], fielddef->m_numModes,modes_offset,&coeffs[m_coeff_offset[eid] + planes_offset*ncoeffs_per_plane]);
}
}
modes_offset += (*m_exp)[0]->GetNumBases();
}
}
}
void Nektar::MultiRegions::ExpListHomogeneous1D::v_FwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
CoeffState  coeffstate 
)
protectedvirtual

Forward transform

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 222 of file ExpListHomogeneous1D.cpp.

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

{
int cnt = 0, cnt1 = 0;
Array<OneD, NekDouble> tmparray;
for(int n = 0; n < m_planes.num_elements(); ++n)
{
m_planes[n]->FwdTrans(inarray+cnt, tmparray = outarray + cnt1,
coeffstate);
cnt += m_planes[n]->GetTotPoints();
cnt1 += m_planes[n]->GetNcoeffs(); // need to skip ncoeffs
}
{
HomogeneousFwdTrans(outarray,outarray,coeffstate);
}
}
void Nektar::MultiRegions::ExpListHomogeneous1D::v_FwdTrans_IterPerExp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Forward transform element by element

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 244 of file ExpListHomogeneous1D.cpp.

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

{
int cnt = 0, cnt1 = 0;
Array<OneD, NekDouble> tmparray;
//spectral element FwdTrans plane by plane
for(int n = 0; n < m_planes.num_elements(); ++n)
{
m_planes[n]->FwdTrans_IterPerExp(inarray+cnt, tmparray = outarray + cnt1);
cnt += m_planes[n]->GetTotPoints();
cnt1 += m_planes[n]->GetNcoeffs();
}
{
HomogeneousFwdTrans(outarray,outarray);
}
}
std::vector< LibUtilities::FieldDefinitionsSharedPtr > Nektar::MultiRegions::ExpListHomogeneous1D::v_GetFieldDefinitions ( void  )
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 572 of file ExpListHomogeneous1D.cpp.

References m_homogeneousBasis, m_lhom, m_planes, and m_transposition.

{
std::vector<LibUtilities::FieldDefinitionsSharedPtr> returnval;
// Set up Homogeneous length details.
Array<OneD,LibUtilities::BasisSharedPtr> HomoBasis(1,m_homogeneousBasis);
std::vector<NekDouble> HomoLen;
HomoLen.push_back(m_lhom);
std::vector<unsigned int> PlanesIDs;
for(int i = 0; i < m_planes.num_elements(); i++)
{
PlanesIDs.push_back(m_transposition->GetPlaneID(i));
}
m_planes[0]->GeneralGetFieldDefinitions(returnval, 1, HomoBasis, HomoLen, PlanesIDs);
return returnval;
}
void Nektar::MultiRegions::ExpListHomogeneous1D::v_GetFieldDefinitions ( std::vector< LibUtilities::FieldDefinitionsSharedPtr > &  fielddef)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 594 of file ExpListHomogeneous1D.cpp.

References m_homogeneousBasis, m_lhom, m_planes, and m_transposition.

{
// Set up Homogeneous length details.
Array<OneD,LibUtilities::BasisSharedPtr> HomoBasis(1,m_homogeneousBasis);
std::vector<NekDouble> HomoLen;
HomoLen.push_back(m_lhom);
std::vector<unsigned int> PlanesIDs;
for(int i = 0; i < m_planes.num_elements(); i++)
{
PlanesIDs.push_back(m_transposition->GetPlaneID(i));
}
// enforce NumHomoDir == 1 by direct call
m_planes[0]->GeneralGetFieldDefinitions(fielddef,1, HomoBasis,HomoLen,PlanesIDs);
}
virtual LibUtilities::BasisSharedPtr Nektar::MultiRegions::ExpListHomogeneous1D::v_GetHomogeneousBasis ( void  )
inlineprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 185 of file ExpListHomogeneous1D.h.

References GetHomogeneousBasis().

{
}
NekDouble Nektar::MultiRegions::ExpListHomogeneous1D::v_GetHomoLen ( void  )
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1068 of file ExpListHomogeneous1D.cpp.

References m_lhom.

{
return m_lhom;
}
virtual int Nektar::MultiRegions::ExpListHomogeneous1D::v_GetNumElmts ( void  )
inlineprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 180 of file ExpListHomogeneous1D.h.

References m_planes.

{
return m_planes[0]->GetExpSize();
}
virtual ExpListSharedPtr& Nektar::MultiRegions::ExpListHomogeneous1D::v_GetPlane ( int  n)
inlineprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 262 of file ExpListHomogeneous1D.h.

References GetPlane().

{
return GetPlane(n);
}
LibUtilities::TranspositionSharedPtr Nektar::MultiRegions::ExpListHomogeneous1D::v_GetTransposition ( void  )
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1063 of file ExpListHomogeneous1D.cpp.

References m_transposition.

{
}
Array< OneD, const unsigned int > Nektar::MultiRegions::ExpListHomogeneous1D::v_GetZIDs ( void  )
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1073 of file ExpListHomogeneous1D.cpp.

References m_transposition.

{
return m_transposition->GetPlanesIDs();
}
void Nektar::MultiRegions::ExpListHomogeneous1D::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 131 of file ExpListHomogeneous1D.cpp.

References Homogeneous1DTrans().

Referenced by HomogeneousBwdTrans().

{
// Backwards trans
Homogeneous1DTrans(inarray,outarray,false,coeffstate,Shuff,UnShuff);
}
void Nektar::MultiRegions::ExpListHomogeneous1D::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 121 of file ExpListHomogeneous1D.cpp.

References Homogeneous1DTrans().

Referenced by HomogeneousFwdTrans().

{
// Forwards trans
Homogeneous1DTrans(inarray,outarray,true,coeffstate,Shuff,UnShuff);
}
void Nektar::MultiRegions::ExpListHomogeneous1D::v_IProductWRTBase ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
CoeffState  coeffstate 
)
protectedvirtual

Inner product

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 307 of file ExpListHomogeneous1D.cpp.

References m_planes.

{
int cnt = 0, cnt1 = 0;
Array<OneD, NekDouble> tmparray;
for(int n = 0; n < m_planes.num_elements(); ++n)
{
m_planes[n]->IProductWRTBase(inarray+cnt, tmparray = outarray + cnt1,coeffstate);
cnt1 += m_planes[n]->GetNcoeffs();
cnt += m_planes[n]->GetTotPoints();
}
}
void Nektar::MultiRegions::ExpListHomogeneous1D::v_IProductWRTBase_IterPerExp ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Inner product element by element

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 324 of file ExpListHomogeneous1D.cpp.

References m_planes.

{
int cnt = 0, cnt1 = 0;
Array<OneD, NekDouble> tmparray;
for(int n = 0; n < m_planes.num_elements(); ++n)
{
m_planes[n]->IProductWRTBase_IterPerExp(inarray+cnt, tmparray = outarray + cnt1);
cnt1 += m_planes[n]->GetNcoeffs();
cnt += m_planes[n]->GetTotPoints();
}
}
void Nektar::MultiRegions::ExpListHomogeneous1D::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 848 of file ExpListHomogeneous1D.cpp.

References ASSERTL0, Nektar::LibUtilities::eFourier, Nektar::LibUtilities::eFourierHalfModeIm, Nektar::LibUtilities::eFourierHalfModeRe, Nektar::LibUtilities::eFourierSingleMode, Nektar::LibUtilities::eXYtoZ, Nektar::LibUtilities::eZtoXY, HomogeneousBwdTrans(), HomogeneousFwdTrans(), Nektar::MultiRegions::ExpList::m_comm, m_homogeneousBasis, m_lhom, m_planes, m_transposition, Nektar::MultiRegions::ExpList::m_WaveSpace, Nektar::NullNekDouble1DArray, sign, and Vmath::Smul().

Referenced by PhysDeriv().

{
int nT_pts = inarray.num_elements(); //number of total points = n. of Fourier points * n. of points per plane (nT_pts)
int nP_pts = nT_pts/m_planes.num_elements(); //number of points per plane = n of Fourier transform required (nP_pts)
Array<OneD, NekDouble> temparray(nT_pts);
Array<OneD, NekDouble> outarray(nT_pts);
Array<OneD, NekDouble> tmp1;
Array<OneD, NekDouble> tmp2;
Array<OneD, NekDouble> tmp3;
for(int i = 0; i < m_planes.num_elements(); i++)
{
m_planes[i]->PhysDeriv(inarray + i*nP_pts ,tmp2 = out_d0 + i*nP_pts , tmp3 = out_d1 + i*nP_pts );
}
if(out_d2 != NullNekDouble1DArray)
{
{
{
temparray = inarray;
}
else
{
HomogeneousFwdTrans(inarray,temparray);
}
NekDouble sign = -1.0;
NekDouble beta;
//Half Mode
{
beta = sign*2*M_PI*(m_transposition->GetK(0))/m_lhom;
Vmath::Smul(nP_pts,beta,temparray,1,outarray,1);
}
{
beta = -sign*2*M_PI*(m_transposition->GetK(0))/m_lhom;
Vmath::Smul(nP_pts,beta,temparray,1,outarray,1);
}
//Fully complex
else
{
for(int i = 0; i < m_planes.num_elements(); i++)
{
beta = -sign*2*M_PI*(m_transposition->GetK(i))/m_lhom;
Vmath::Smul(nP_pts,beta,tmp1 = temparray + i*nP_pts,1,tmp2 = outarray + (i-int(sign))*nP_pts,1);
sign = -1.0*sign;
}
}
{
out_d2 = outarray;
}
else
{
HomogeneousBwdTrans(outarray,out_d2);
}
}
else
{
ASSERTL0(m_comm->GetColumnComm()->GetSize() == 1,"Parallelisation in the homogeneous direction implemented just for Fourier basis");
{
ASSERTL0(false,"Semi-phyisical time-stepping not implemented yet for non-Fourier basis");
}
else
{
StdRegions::StdSegExp StdSeg(m_homogeneousBasis->GetBasisKey());
m_transposition->Transpose(inarray,temparray,false,LibUtilities::eXYtoZ);
for(int i = 0; i < nP_pts; i++)
{
StdSeg.PhysDeriv(temparray + i*m_planes.num_elements(), tmp2 = outarray + i*m_planes.num_elements());
}
m_transposition->Transpose(outarray,out_d2,false,LibUtilities::eZtoXY);
Vmath::Smul(nT_pts,2.0/m_lhom,out_d2,1,out_d2,1);
}
}
}
}
void Nektar::MultiRegions::ExpListHomogeneous1D::v_PhysDeriv ( Direction  edir,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  out_d 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 948 of file ExpListHomogeneous1D.cpp.

References ASSERTL0, Nektar::LibUtilities::eFourier, Nektar::LibUtilities::eFourierHalfModeIm, Nektar::LibUtilities::eFourierHalfModeRe, Nektar::LibUtilities::eFourierSingleMode, Nektar::LibUtilities::eXYtoZ, Nektar::LibUtilities::eZtoXY, HomogeneousBwdTrans(), HomogeneousFwdTrans(), Nektar::MultiRegions::ExpList::m_comm, m_homogeneousBasis, m_lhom, m_planes, m_transposition, Nektar::MultiRegions::ExpList::m_WaveSpace, sign, and Vmath::Smul().

{
int nT_pts = inarray.num_elements(); //number of total points = n. of Fourier points * n. of points per plane (nT_pts)
int nP_pts = nT_pts/m_planes.num_elements(); //number of points per plane = n of Fourier transform required (nP_pts)
int dir= (int)edir;
Array<OneD, NekDouble> temparray(nT_pts);
Array<OneD, NekDouble> outarray(nT_pts);
Array<OneD, NekDouble> tmp1;
Array<OneD, NekDouble> tmp2;
if (dir < 2)
{
for(int i=0; i < m_planes.num_elements(); i++)
{
m_planes[i]->PhysDeriv(edir, inarray + i*nP_pts ,tmp2 = out_d + i*nP_pts);
}
}
else
{
{
{
temparray = inarray;
}
else
{
HomogeneousFwdTrans(inarray,temparray);
}
NekDouble sign = -1.0;
NekDouble beta;
//HalfMode
{
beta = 2*M_PI*(m_transposition->GetK(0))/m_lhom;
Vmath::Smul(nP_pts,beta,temparray,1,outarray,1);
}
{
beta = -2*M_PI*(m_transposition->GetK(0))/m_lhom;
Vmath::Smul(nP_pts,beta,temparray,1,outarray,1);
}
//Fully complex
else
{
for(int i = 0; i < m_planes.num_elements(); i++)
{
beta = -sign*2*M_PI*(m_transposition->GetK(i))/m_lhom;
Vmath::Smul(nP_pts,beta,tmp1 = temparray + i*nP_pts,1,tmp2 = outarray + (i-int(sign))*nP_pts,1);
sign = -1.0*sign;
}
}
{
out_d = outarray;
}
else
{
HomogeneousBwdTrans(outarray,out_d);
}
}
else
{
ASSERTL0(m_comm->GetColumnComm()->GetSize() == 1,"Parallelisation in the homogeneous direction implemented just for Fourier basis");
{
ASSERTL0(false,"Semi-phyisical time-stepping not implemented yet for non-Fourier basis");
}
else
{
StdRegions::StdSegExp StdSeg(m_homogeneousBasis->GetBasisKey());
m_transposition->Transpose(inarray,temparray,false,LibUtilities::eXYtoZ);
for(int i = 0; i < nP_pts; i++)
{
StdSeg.PhysDeriv(temparray + i*m_planes.num_elements(), tmp2 = outarray + i*m_planes.num_elements());
}
m_transposition->Transpose(outarray,out_d,false,LibUtilities::eZtoXY);
Vmath::Smul(nT_pts,2.0/m_lhom,out_d,1,out_d,1);
}
}
}
}
void Nektar::MultiRegions::ExpListHomogeneous1D::v_PhysGalerkinProjection1DScaled ( const NekDouble  scale,
const Array< OneD, NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 831 of file ExpListHomogeneous1D.cpp.

References ASSERTL1, and m_planes.

{
int cnt,cnt1;
Array<OneD, NekDouble> tmparray;
cnt = m_planes[0]->Get1DScaledTotPoints(scale);
cnt1 = m_planes[0]->GetTotPoints();
ASSERTL1(m_planes.num_elements()*cnt <= inarray.num_elements(),"size of outarray does not match internal estimage");
for(int i = 0; i < m_planes.num_elements(); i++)
{
m_planes[i]->PhysGalerkinProjection1DScaled(scale,inarray+i*cnt,
tmparray = outarray+i*cnt1);
}
}
void Nektar::MultiRegions::ExpListHomogeneous1D::v_PhysInterp1DScaled ( const NekDouble  scale,
const Array< OneD, NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 812 of file ExpListHomogeneous1D.cpp.

References ASSERTL1, and m_planes.

{
int cnt,cnt1;
Array<OneD, NekDouble> tmparray;
cnt = m_planes[0]->GetTotPoints();
cnt1 = m_planes[0]->Get1DScaledTotPoints(scale);
ASSERTL1(m_planes.num_elements()*cnt1 <= outarray.num_elements(),"size of outarray does not match internal estimage");
for(int i = 0; i < m_planes.num_elements(); i++)
{
m_planes[i]->PhysInterp1DScaled(scale,inarray+i*cnt,
tmparray = outarray+i*cnt1);
}
}
virtual void Nektar::MultiRegions::ExpListHomogeneous1D::v_SetHomo1DSpecVanVisc ( Array< OneD, NekDouble visc)
inlineprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 175 of file ExpListHomogeneous1D.h.

References m_specVanVisc.

{
m_specVanVisc = visc;
}
void Nektar::MultiRegions::ExpListHomogeneous1D::v_WriteVtkPieceData ( std::ofstream &  outfile,
int  expansion,
std::string  var 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 789 of file ExpListHomogeneous1D.cpp.

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

{
int i;
int nq = (*m_exp)[expansion]->GetTotPoints();
int npoints_per_plane = m_planes[0]->GetTotPoints();
// printing the fields of that zone
outfile << " <DataArray type=\"Float32\" Name=\""
<< var << "\">" << endl;
outfile << " ";
for (int n = 0; n < m_planes.num_elements(); ++n)
{
const Array<OneD, NekDouble> phys = m_phys + m_phys_offset[expansion] + n*npoints_per_plane;
for(i = 0; i < nq; ++i)
{
outfile << (fabs(phys[i]) < NekConstants::kNekZeroTol ? 0 : phys[i]) << " ";
}
}
outfile << endl;
outfile << " </DataArray>" << endl;
}

Member Data Documentation

bool Nektar::MultiRegions::ExpListHomogeneous1D::m_dealiasing
private

Definition at line 272 of file ExpListHomogeneous1D.h.

Referenced by ExpListHomogeneous1D().

LibUtilities::NektarFFTSharedPtr Nektar::MultiRegions::ExpListHomogeneous1D::m_FFT
protected

Definition at line 142 of file ExpListHomogeneous1D.h.

Referenced by ExpListHomogeneous1D(), and Homogeneous1DTrans().

LibUtilities::NektarFFTSharedPtr Nektar::MultiRegions::ExpListHomogeneous1D::m_FFT_deal
protected

Definition at line 144 of file ExpListHomogeneous1D.h.

Referenced by ExpListHomogeneous1D(), and v_DealiasedProd().

Homo1DBlockMatrixMapShPtr Nektar::MultiRegions::ExpListHomogeneous1D::m_homogeneous1DBlockMat
protected

Definition at line 154 of file ExpListHomogeneous1D.h.

Referenced by GetHomogeneous1DBlockMatrix().

LibUtilities::BasisSharedPtr Nektar::MultiRegions::ExpListHomogeneous1D::m_homogeneousBasis
protected

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

Definition at line 152 of file ExpListHomogeneous1D.h.

Referenced by Nektar::MultiRegions::DisContField3DHomogeneous1D::EvaluateBoundaryConditions(), ExpListHomogeneous1D(), Nektar::MultiRegions::ExpList3DHomogeneous1D::GenExpList3DHomogeneous1D(), GenHomogeneous1DBlockMatrix(), Nektar::MultiRegions::ExpList2DHomogeneous1D::GetCoords(), Nektar::MultiRegions::ExpList3DHomogeneous1D::GetCoords(), GetHomogeneousBasis(), Homogeneous1DTrans(), v_DealiasedProd(), v_ExtractDataToCoeffs(), Nektar::MultiRegions::ExpList2DHomogeneous1D::v_GetCoords(), Nektar::MultiRegions::ExpList3DHomogeneous1D::v_GetCoords(), v_GetFieldDefinitions(), Nektar::MultiRegions::ContField3DHomogeneous1D::v_HelmSolve(), Nektar::MultiRegions::ExpList3DHomogeneous1D::v_L2(), and v_PhysDeriv().

NekDouble Nektar::MultiRegions::ExpListHomogeneous1D::m_lhom
protected

Width of homogeneous direction.

Definition at line 153 of file ExpListHomogeneous1D.h.

Referenced by Nektar::MultiRegions::DisContField3DHomogeneous1D::EvaluateBoundaryConditions(), Nektar::MultiRegions::ExpList2DHomogeneous1D::GetCoords(), Nektar::MultiRegions::ExpList3DHomogeneous1D::GetCoords(), Nektar::MultiRegions::ExpList2DHomogeneous1D::v_GetCoords(), Nektar::MultiRegions::ExpList3DHomogeneous1D::v_GetCoords(), v_GetFieldDefinitions(), v_GetHomoLen(), Nektar::MultiRegions::ContField3DHomogeneous1D::v_HelmSolve(), Nektar::MultiRegions::DisContField3DHomogeneous1D::v_HelmSolve(), Nektar::MultiRegions::ExpList3DHomogeneous1D::v_L2(), and v_PhysDeriv().

int Nektar::MultiRegions::ExpListHomogeneous1D::m_padsize
private

Definition at line 273 of file ExpListHomogeneous1D.h.

Referenced by ExpListHomogeneous1D(), and v_DealiasedProd().

Array<OneD, ExpListSharedPtr> Nektar::MultiRegions::ExpListHomogeneous1D::m_planes
protected

Definition at line 155 of file ExpListHomogeneous1D.h.

Referenced by Nektar::MultiRegions::ContField3DHomogeneous1D::ContField3DHomogeneous1D(), Nektar::MultiRegions::DisContField3DHomogeneous1D::DisContField3DHomogeneous1D(), Nektar::MultiRegions::DisContField3DHomogeneous1D::EvaluateBoundaryConditions(), Nektar::MultiRegions::ExpList2DHomogeneous1D::ExpList2DHomogeneous1D(), Nektar::MultiRegions::ExpList3DHomogeneous1D::ExpList3DHomogeneous1D(), ExpListHomogeneous1D(), Nektar::MultiRegions::ExpList3DHomogeneous1D::GenExpList3DHomogeneous1D(), GenHomogeneous1DBlockMatrix(), Nektar::MultiRegions::DisContField3DHomogeneous1D::GetBCValues(), Nektar::MultiRegions::DisContField3DHomogeneous1D::GetBoundaryToElmtMap(), Nektar::MultiRegions::ExpList2DHomogeneous1D::GetCoords(), Nektar::MultiRegions::ExpList3DHomogeneous1D::GetCoords(), GetPlane(), Homogeneous1DTrans(), Nektar::MultiRegions::DisContField3DHomogeneous1D::NormVectorIProductWRTBase(), Nektar::MultiRegions::ExpList2DHomogeneous1D::SetCoeffPhys(), Nektar::MultiRegions::ExpList3DHomogeneous1D::SetCoeffPhys(), Nektar::MultiRegions::DisContField3DHomogeneous1D::SetupBoundaryConditions(), Nektar::MultiRegions::DisContField3DHomogeneous1D::SetUpDG(), v_AppendFieldData(), v_BwdTrans(), v_BwdTrans_IterPerExp(), v_DealiasedProd(), v_ExtractCoeffsToCoeffs(), v_ExtractDataToCoeffs(), Nektar::MultiRegions::DisContField3DHomogeneous1D::v_ExtractTracePhys(), v_FwdTrans(), v_FwdTrans_IterPerExp(), Nektar::MultiRegions::ExpList2DHomogeneous1D::v_GetCoords(), Nektar::MultiRegions::ExpList3DHomogeneous1D::v_GetCoords(), v_GetFieldDefinitions(), Nektar::MultiRegions::ExpList2DHomogeneous1D::v_GetNormals(), v_GetNumElmts(), Nektar::MultiRegions::ExpList3DHomogeneous1D::v_GetPeriodicEntities(), Nektar::MultiRegions::DisContField3DHomogeneous1D::v_GetTraceMap(), Nektar::MultiRegions::ContField3DHomogeneous1D::v_GlobalToLocal(), Nektar::MultiRegions::ContField3DHomogeneous1D::v_HelmSolve(), Nektar::MultiRegions::DisContField3DHomogeneous1D::v_HelmSolve(), Nektar::MultiRegions::ExpList3DHomogeneous1D::v_HomogeneousEnergy(), Nektar::MultiRegions::ContField3DHomogeneous1D::v_ImposeDirichletConditions(), v_IProductWRTBase(), v_IProductWRTBase_IterPerExp(), Nektar::MultiRegions::ExpList3DHomogeneous1D::v_L2(), Nektar::MultiRegions::ContField3DHomogeneous1D::v_LocalToGlobal(), v_PhysDeriv(), v_PhysGalerkinProjection1DScaled(), v_PhysInterp1DScaled(), Nektar::MultiRegions::ContField3DHomogeneous1D::v_SmoothField(), Nektar::MultiRegions::ExpList3DHomogeneous1D::v_WriteTecplotConnectivity(), Nektar::MultiRegions::ExpList2DHomogeneous1D::v_WriteTecplotZone(), v_WriteVtkPieceData(), Nektar::MultiRegions::ExpList3DHomogeneous1D::v_WriteVtkPieceHeader(), and Nektar::MultiRegions::ExpList2DHomogeneous1D::v_WriteVtkPieceHeader().

Array<OneD, NekDouble> Nektar::MultiRegions::ExpListHomogeneous1D::m_specVanVisc
private

Spectral vanishing Viscosity coefficient for stabilisation.

Definition at line 276 of file ExpListHomogeneous1D.h.

Referenced by GetSpecVanVisc(), and v_SetHomo1DSpecVanVisc().

Array<OneD,NekDouble> Nektar::MultiRegions::ExpListHomogeneous1D::m_tmpIN
protected

Definition at line 146 of file ExpListHomogeneous1D.h.

Referenced by Homogeneous1DTrans().

Array<OneD,NekDouble> Nektar::MultiRegions::ExpListHomogeneous1D::m_tmpOUT
protected

Definition at line 147 of file ExpListHomogeneous1D.h.

Referenced by Homogeneous1DTrans().

LibUtilities::TranspositionSharedPtr Nektar::MultiRegions::ExpListHomogeneous1D::m_transposition

Definition at line 136 of file ExpListHomogeneous1D.h.

Referenced by Nektar::MultiRegions::ContField3DHomogeneous1D::ContField3DHomogeneous1D(), Nektar::MultiRegions::DisContField3DHomogeneous1D::EvaluateBoundaryConditions(), ExpListHomogeneous1D(), Nektar::MultiRegions::ExpList2DHomogeneous1D::GetCoords(), Nektar::MultiRegions::ExpList3DHomogeneous1D::GetCoords(), Homogeneous1DTrans(), v_DealiasedProd(), v_ExtractDataToCoeffs(), Nektar::MultiRegions::ExpList2DHomogeneous1D::v_GetCoords(), Nektar::MultiRegions::ExpList3DHomogeneous1D::v_GetCoords(), v_GetFieldDefinitions(), v_GetTransposition(), v_GetZIDs(), Nektar::MultiRegions::ContField3DHomogeneous1D::v_HelmSolve(), Nektar::MultiRegions::DisContField3DHomogeneous1D::v_HelmSolve(), Nektar::MultiRegions::ExpList3DHomogeneous1D::v_L2(), and v_PhysDeriv().

bool Nektar::MultiRegions::ExpListHomogeneous1D::m_useFFT
protected

FFT variables.

Definition at line 141 of file ExpListHomogeneous1D.h.

Referenced by ExpListHomogeneous1D(), Homogeneous1DTrans(), and Nektar::MultiRegions::DisContField3DHomogeneous1D::SetupBoundaryConditions().