Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Protected Member Functions | Private Member Functions | List of all members
Nektar::MultiRegions::ExpList3DHomogeneous1D Class Reference

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

#include <ExpList3DHomogeneous1D.h>

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

Public Member Functions

 ExpList3DHomogeneous1D ()
 Default constructor. More...
 
 ExpList3DHomogeneous1D (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::BasisKey &HomoBasis, const NekDouble lhom, const bool useFFT, const bool dealiasing)
 
 ExpList3DHomogeneous1D (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::BasisKey &HomoBasis, const NekDouble lhom, const bool useFFT, const bool dealiasing, const SpatialDomains::MeshGraphSharedPtr &graph2D, const std::string &var="DefaultVar")
 Sets up a list of local expansions based on an input mesh. More...
 
 ExpList3DHomogeneous1D (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::BasisKey &HomoBasis, const NekDouble lhom, const bool useFFT, const bool dealiasing, const SpatialDomains::ExpansionMap &expansions)
 Sets up a list of local expansions based on an mesh expansion. More...
 
 ExpList3DHomogeneous1D (const ExpList3DHomogeneous1D &In, const bool DeclarePlanesSetCoeffPhys=true)
 Copy constructor. More...
 
 ExpList3DHomogeneous1D (const ExpList3DHomogeneous1D &In, const std::vector< unsigned int > &eIDs, const bool DeclarePlanesSetCoeffPhys=true)
 Constructor copying only elements defined in eIds. More...
 
virtual ~ExpList3DHomogeneous1D ()
 Destructor. More...
 
void GetCoords (Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1=NullNekDouble1DArray, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
 This function calculates the coordinates of all the elemental quadrature points $\boldsymbol{x}_i$. More...
 
void GetCoords (const int eid, Array< OneD, NekDouble > &xc0, Array< OneD, NekDouble > &xc1, Array< OneD, NekDouble > &xc2)
 
- Public Member Functions inherited from Nektar::MultiRegions::ExpListHomogeneous1D
 ExpListHomogeneous1D ()
 Default constructor. More...
 
 ExpListHomogeneous1D (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::BasisKey &HomoBasis, const NekDouble lz, const bool useFFT, const bool dealiasing)
 
 ExpListHomogeneous1D (const ExpListHomogeneous1D &In)
 Copy constructor. More...
 
 ExpListHomogeneous1D (const ExpListHomogeneous1D &In, const std::vector< unsigned int > &eIDs)
 
virtual ~ExpListHomogeneous1D ()
 Destructor. More...
 
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)
 
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. 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)
 

Protected Member Functions

void SetCoeffPhys (void)
 Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m_coeff and m_phys. More...
 
virtual void v_GetCoords (Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1, Array< OneD, NekDouble > &coord_2)
 
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 > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
 
virtual Array< OneD, const
NekDouble
v_HomogeneousEnergy (void)
 
virtual void v_GetPeriodicEntities (PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces)
 
- Protected Member Functions inherited from Nektar::MultiRegions::ExpListHomogeneous1D
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. More...
 
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::ostream &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. 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_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_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 NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray)
 
virtual Array< OneD, const
unsigned int > 
v_GetYIDs (void)
 
virtual void v_ClearGlobalLinSysManager (void)
 
void ExtractFileBCs (const std::string &fileName, const std::string &varName, const boost::shared_ptr< ExpList > locExpList)
 

Private Member Functions

void GenExpList3DHomogeneous1D (const SpatialDomains::ExpansionMap &expansions)
 

Additional Inherited Members

- Public Attributes inherited from Nektar::MultiRegions::ExpListHomogeneous1D
LibUtilities::TranspositionSharedPtr m_transposition
 
LibUtilities::CommSharedPtr m_StripZcomm
 
- Public Attributes inherited from Nektar::MultiRegions::ExpList
ExpansionType m_expType
 
- Static Protected Member Functions inherited from Nektar::MultiRegions::ExpList
static
SpatialDomains::BoundaryConditionShPtr 
GetBoundaryCondition (const SpatialDomains::BoundaryConditionCollection &collection, unsigned int index, const std::string &variable)
 
- Protected Attributes inherited from Nektar::MultiRegions::ExpListHomogeneous1D
bool m_useFFT
 FFT variables. More...
 
LibUtilities::NektarFFTSharedPtr m_FFT
 
LibUtilities::NektarFFTSharedPtr m_FFT_deal
 
Array< OneD, NekDoublem_tmpIN
 
Array< OneD, NekDoublem_tmpOUT
 
LibUtilities::BasisSharedPtr m_homogeneousBasis
 Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m_coeff and m_phys. More...
 
NekDouble m_lhom
 Width of homogeneous direction. More...
 
Homo1DBlockMatrixMapShPtr m_homogeneous1DBlockMat
 
Array< OneD, ExpListSharedPtrm_planes
 
boost::unordered_map< int, int > m_zIdToPlane
 
- 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...
 

Detailed Description

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

Definition at line 61 of file ExpList3DHomogeneous1D.h.

Constructor & Destructor Documentation

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

Default constructor.

Definition at line 47 of file ExpList3DHomogeneous1D.cpp.

References Nektar::MultiRegions::e3DH1D, and Nektar::MultiRegions::ExpList::SetExpType().

47  :
49  {
51  }
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:253
Nektar::MultiRegions::ExpList3DHomogeneous1D::ExpList3DHomogeneous1D ( const LibUtilities::SessionReaderSharedPtr pSession,
const LibUtilities::BasisKey HomoBasis,
const NekDouble  lhom,
const bool  useFFT,
const bool  dealiasing 
)

Definition at line 53 of file ExpList3DHomogeneous1D.cpp.

References Nektar::MultiRegions::e3DH1D, and Nektar::MultiRegions::ExpList::SetExpType().

58  :
59  ExpListHomogeneous1D(pSession,HomoBasis,lhom,useFFT,dealiasing)
60  {
62  }
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:253
Nektar::MultiRegions::ExpList3DHomogeneous1D::ExpList3DHomogeneous1D ( const LibUtilities::SessionReaderSharedPtr pSession,
const LibUtilities::BasisKey HomoBasis,
const NekDouble  lhom,
const bool  useFFT,
const bool  dealiasing,
const SpatialDomains::MeshGraphSharedPtr graph2D,
const std::string &  var = "DefaultVar" 
)

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

Definition at line 65 of file ExpList3DHomogeneous1D.cpp.

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

72  :
73  ExpListHomogeneous1D(pSession,HomoBasis,lhom,useFFT,dealiasing)
74  {
76 
77  GenExpList3DHomogeneous1D(graph2D->GetExpansions(var));
78  }
void GenExpList3DHomogeneous1D(const SpatialDomains::ExpansionMap &expansions)
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:253
Nektar::MultiRegions::ExpList3DHomogeneous1D::ExpList3DHomogeneous1D ( const LibUtilities::SessionReaderSharedPtr pSession,
const LibUtilities::BasisKey HomoBasis,
const NekDouble  lhom,
const bool  useFFT,
const bool  dealiasing,
const SpatialDomains::ExpansionMap expansions 
)

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

Definition at line 81 of file ExpList3DHomogeneous1D.cpp.

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

87  :
88  ExpListHomogeneous1D(pSession,HomoBasis,lhom,useFFT,dealiasing)
89  {
91 
92  GenExpList3DHomogeneous1D(expansions);
93  }
void GenExpList3DHomogeneous1D(const SpatialDomains::ExpansionMap &expansions)
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:253
Nektar::MultiRegions::ExpList3DHomogeneous1D::ExpList3DHomogeneous1D ( const ExpList3DHomogeneous1D In,
const bool  DeclarePlanesSetCoeffPhys = true 
)

Copy constructor.

Parameters
InExpList3DHomogeneous1D object to copy.

Definition at line 133 of file ExpList3DHomogeneous1D.cpp.

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

133  :
135  {
137 
138  if(DeclarePlanesSetCoeffPhys)
139  {
140  bool False = false;
141  ExpList2DSharedPtr zero_plane = boost::dynamic_pointer_cast<ExpList2D> (In.m_planes[0]);
142 
143  for(int n = 0; n < m_planes.num_elements(); ++n)
144  {
146  }
147 
148  SetCoeffPhys();
149  }
150  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
void SetCoeffPhys(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
Array< OneD, ExpListSharedPtr > m_planes
boost::shared_ptr< ExpList2D > ExpList2DSharedPtr
Shared pointer to an ExpList2D object.
Definition: ExpList2D.h:49
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:253
Nektar::MultiRegions::ExpList3DHomogeneous1D::ExpList3DHomogeneous1D ( const ExpList3DHomogeneous1D In,
const std::vector< unsigned int > &  eIDs,
const bool  DeclarePlanesSetCoeffPhys = true 
)

Constructor copying only elements defined in eIds.

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

Definition at line 156 of file ExpList3DHomogeneous1D.cpp.

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

158  :
159  ExpListHomogeneous1D(In, eIDs)
160  {
162 
163  if(DeclarePlanesSetCoeffPhys)
164  {
165  bool False = false;
166  std::vector<unsigned int> eIDsPlane;
167  int nel = eIDs.size()/m_planes.num_elements();
168 
169  for (int i = 0; i < nel; ++i)
170  {
171  eIDsPlane.push_back(eIDs[i]);
172  }
173 
174  ExpList2DSharedPtr zero_plane_old =
175  boost::dynamic_pointer_cast<ExpList2D> (In.m_planes[0]);
176 
177  ExpList2DSharedPtr zero_plane =
178  MemoryManager<ExpList2D>::AllocateSharedPtr(*(zero_plane_old), eIDsPlane);
179 
180  for(int n = 0; n < m_planes.num_elements(); ++n)
181  {
183  }
184 
185  SetCoeffPhys();
186  }
187  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
void SetCoeffPhys(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
Array< OneD, ExpListSharedPtr > m_planes
boost::shared_ptr< ExpList2D > ExpList2DSharedPtr
Shared pointer to an ExpList2D object.
Definition: ExpList2D.h:49
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:253
Nektar::MultiRegions::ExpList3DHomogeneous1D::~ExpList3DHomogeneous1D ( )
virtual

Destructor.

Destructor

Definition at line 192 of file ExpList3DHomogeneous1D.cpp.

193  {
194  }

Member Function Documentation

void Nektar::MultiRegions::ExpList3DHomogeneous1D::GenExpList3DHomogeneous1D ( const SpatialDomains::ExpansionMap expansions)
private

Definition at line 95 of file ExpList3DHomogeneous1D.cpp.

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

Referenced by ExpList3DHomogeneous1D().

96  {
97  int n,j,nel;
98  bool False = false;
99  ExpList2DSharedPtr plane_zero;
100 
101  // note that nzplanes can be larger than nzmodes
102  m_planes[0] = plane_zero = MemoryManager<ExpList2D>::AllocateSharedPtr(m_session, expansions, False);
103 
105  nel = m_planes[0]->GetExpSize();
106 
107  for(j = 0; j < nel; ++j)
108  {
109  (*m_exp).push_back(m_planes[0]->GetExp(j));
110  }
111 
112  for(n = 1; n < m_planes.num_elements(); ++n)
113  {
115  for(j = 0; j < nel; ++j)
116  {
117  (*m_exp).push_back((*m_exp)[j]);
118  }
119  }
120 
121  // Setup Default optimisation information.
122  nel = GetExpSize();
125 
126  SetCoeffPhys();
127  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:1001
const boost::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:1920
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1899
void SetCoeffPhys(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:977
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:910
Array< OneD, ExpListSharedPtr > m_planes
boost::shared_ptr< ExpList2D > ExpList2DSharedPtr
Shared pointer to an ExpList2D object.
Definition: ExpList2D.h:49
void Nektar::MultiRegions::ExpList3DHomogeneous1D::GetCoords ( Array< OneD, NekDouble > &  coord_0,
Array< OneD, NekDouble > &  coord_1 = NullNekDouble1DArray,
Array< OneD, NekDouble > &  coord_2 = NullNekDouble1DArray 
)
inline

This function calculates the coordinates of all the elemental quadrature points $\boldsymbol{x}_i$.

Definition at line 158 of file ExpList3DHomogeneous1D.h.

References v_GetCoords().

Referenced by v_WriteVtkPieceHeader().

162  {
163  v_GetCoords(coord_0,coord_1,coord_2);
164  }
virtual void v_GetCoords(Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1, Array< OneD, NekDouble > &coord_2)
void Nektar::MultiRegions::ExpList3DHomogeneous1D::GetCoords ( const int  eid,
Array< OneD, NekDouble > &  xc0,
Array< OneD, NekDouble > &  xc1,
Array< OneD, NekDouble > &  xc2 
)

Definition at line 231 of file ExpList3DHomogeneous1D.cpp.

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

235  {
236  int n;
237  Array<OneD, NekDouble> tmp_xc;
238  int nzplanes = m_planes.num_elements();
239  int npoints = GetTotPoints(eid);
240 
241  (*m_exp)[eid]->GetCoords(xc0,xc1);
242 
243  // Fill z-direction
244  Array<OneD, const NekDouble> pts = m_homogeneousBasis->GetZ();
245  Array<OneD, NekDouble> local_pts(m_planes.num_elements());
246 
247  for(n = 0; n < m_planes.num_elements(); n++)
248  {
249  local_pts[n] = pts[m_transposition->GetPlaneID(n)];
250  }
251 
252  Array<OneD, NekDouble> z(nzplanes);
253 
254  Vmath::Smul(nzplanes,m_lhom/2.0,local_pts,1,z,1);
255  Vmath::Sadd(nzplanes,m_lhom/2.0,z,1,z,1);
256 
257  for(n = 0; n < nzplanes; ++n)
258  {
259  Vmath::Fill(npoints,z[n],tmp_xc = xc2 + npoints*n,1);
260  if(n)
261  {
262  Vmath::Vcopy(npoints,xc0,1,tmp_xc = xc0+npoints*n,1);
263  Vmath::Vcopy(npoints,xc1,1,tmp_xc = xc1+npoints*n,1);
264  }
265  }
266  }
LibUtilities::TranspositionSharedPtr m_transposition
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:46
NekDouble m_lhom
Width of homogeneous direction.
LibUtilities::BasisSharedPtr m_homogeneousBasis
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
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
int GetTotPoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1456
Array< OneD, ExpListSharedPtr > m_planes
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition: Vmath.cpp:301
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
void Nektar::MultiRegions::ExpList3DHomogeneous1D::SetCoeffPhys ( void  )
protected

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

Definition at line 196 of file ExpList3DHomogeneous1D.cpp.

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

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

197  {
198  int i,n,cnt;
199  int ncoeffs_per_plane = m_planes[0]->GetNcoeffs();
200  int npoints_per_plane = m_planes[0]->GetTotPoints();
201 
202  int nzplanes = m_planes.num_elements();
203 
204  // Set total coefficients and points
205  m_ncoeffs = ncoeffs_per_plane*nzplanes;
206  m_npoints = npoints_per_plane*nzplanes;
207 
208  m_coeffs = Array<OneD, NekDouble> (m_ncoeffs,0.0);
209  m_phys = Array<OneD, NekDouble> (m_npoints,0.0);
210 
211  int nel = m_planes[0]->GetExpSize();
212  m_coeff_offset = Array<OneD,int>(nel*nzplanes);
213  m_phys_offset = Array<OneD,int>(nel*nzplanes);
214  m_offset_elmt_id = Array<OneD,int>(nel*nzplanes);
215  Array<OneD, NekDouble> tmparray;
216 
217  for(cnt = n = 0; n < nzplanes; ++n)
218  {
219  m_planes[n]->SetCoeffsArray(tmparray= m_coeffs + ncoeffs_per_plane*n);
220  m_planes[n]->SetPhysArray(tmparray = m_phys + npoints_per_plane*n);
221 
222  for(i = 0; i < nel; ++i)
223  {
224  m_coeff_offset[cnt] = m_planes[n]->GetCoeff_Offset(i) + n*ncoeffs_per_plane;
225  m_phys_offset[cnt] = m_planes[n]->GetPhys_Offset(i) + n*npoints_per_plane;
226  m_offset_elmt_id[cnt++] = m_planes[n]->GetOffset_Elmt_Id(i) + n*nel;
227  }
228  }
229  }
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:956
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:939
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:988
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:991
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:917
Array< OneD, int > m_offset_elmt_id
Array containing the element id m_offset_elmt_id[n] that the n^th consecutive block of data in m_coef...
Definition: ExpList.h:999
Array< OneD, ExpListSharedPtr > m_planes
void Nektar::MultiRegions::ExpList3DHomogeneous1D::v_GetCoords ( Array< OneD, NekDouble > &  xc0,
Array< OneD, NekDouble > &  xc1,
Array< OneD, NekDouble > &  xc2 
)
protectedvirtual

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

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

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 284 of file ExpList3DHomogeneous1D.cpp.

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

Referenced by GetCoords().

287  {
288  int n;
289  Array<OneD, NekDouble> tmp_xc;
290  int nzplanes = m_planes.num_elements();
291  int npoints = m_planes[0]->GetTotPoints();
292 
293  m_planes[0]->GetCoords(xc0,xc1);
294 
295  // Fill z-direction
296  Array<OneD, const NekDouble> pts = m_homogeneousBasis->GetZ();
297 
298  Array<OneD, NekDouble> local_pts(m_planes.num_elements());
299 
300  for(n = 0; n < m_planes.num_elements(); n++)
301  {
302  local_pts[n] = pts[m_transposition->GetPlaneID(n)];
303  }
304 
305  Array<OneD, NekDouble> z(nzplanes);
306 
307  Vmath::Smul(nzplanes,m_lhom/2.0,local_pts,1,z,1);
308  Vmath::Sadd(nzplanes,m_lhom/2.0,z,1,z,1);
309 
310  for(n = 0; n < nzplanes; ++n)
311  {
312  Vmath::Fill(npoints,z[n],tmp_xc = xc2 + npoints*n,1);
313  if(n)
314  {
315  Vmath::Vcopy(npoints,xc0,1,tmp_xc = xc0+npoints*n,1);
316  Vmath::Vcopy(npoints,xc1,1,tmp_xc = xc1+npoints*n,1);
317  }
318  }
319  }
LibUtilities::TranspositionSharedPtr m_transposition
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:46
NekDouble m_lhom
Width of homogeneous direction.
LibUtilities::BasisSharedPtr m_homogeneousBasis
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
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
Array< OneD, ExpListSharedPtr > m_planes
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition: Vmath.cpp:301
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
virtual void Nektar::MultiRegions::ExpList3DHomogeneous1D::v_GetPeriodicEntities ( PeriodicMap periodicVerts,
PeriodicMap periodicEdges,
PeriodicMap periodicFaces 
)
inlineprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 142 of file ExpList3DHomogeneous1D.h.

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

147  {
148  m_planes[0]->GetPeriodicEntities(periodicVerts,periodicEdges);
149  }
Array< OneD, ExpListSharedPtr > m_planes
Array< OneD, const NekDouble > Nektar::MultiRegions::ExpList3DHomogeneous1D::v_HomogeneousEnergy ( void  )
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 515 of file ExpList3DHomogeneous1D.cpp.

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

516  {
517  Array<OneD, NekDouble> energy(m_planes.num_elements()/2);
518  NekDouble area = 0.0;
519 
520  // Calculate total area of elements.
521  for (int n = 0; n < m_planes[0]->GetExpSize(); ++n)
522  {
523  Array<OneD, NekDouble> inarray(m_planes[0]->GetExp(n)->GetTotPoints(), 1.0);
524  area += m_planes[0]->GetExp(n)->Integral(inarray);
525  }
526 
527  m_comm->GetRowComm()->AllReduce(area, LibUtilities::ReduceSum);
528 
529  // Calculate L2 norm of real/imaginary planes.
530  for (int n = 0; n < m_planes.num_elements(); n += 2)
531  {
532  NekDouble err;
533 
534  energy[n/2] = 0;
535 
536  for(int i = 0; i < m_planes[n]->GetExpSize(); ++i)
537  {
538  StdRegions::StdExpansionSharedPtr exp = m_planes[n]->GetExp(i);
539  Array<OneD, NekDouble> phys(exp->GetTotPoints());
540  exp->BwdTrans(m_planes[n]->GetCoeffs()+m_planes[n]->GetCoeff_Offset(i),
541  phys);
542  err = exp->L2(phys);
543  energy[n/2] += err*err;
544 
545  exp = m_planes[n+1]->GetExp(i);
546  exp->BwdTrans(m_planes[n+1]->GetCoeffs()+m_planes[n+1]->GetCoeff_Offset(i),
547  phys);
548  err = exp->L2(phys);
549  energy[n/2] += err*err;
550  }
551 
552  m_comm->GetRowComm()->AllReduce(energy[n/2], LibUtilities::ReduceSum);
553  energy[n/2] /= 2.0*area;
554  }
555 
556  return energy;
557  }
const Array< OneD, const NekDouble > & GetCoeffs() const
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expa...
Definition: ExpList.h:1837
int GetCoeff_Offset(int n) const
Get the start offset position for a global list of m_coeffs correspoinding to element n...
Definition: ExpList.h:1929
const boost::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:1920
int GetTotPoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1456
Array< OneD, ExpListSharedPtr > m_planes
double NekDouble
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: ExpList.h:907
boost::shared_ptr< StdExpansion > StdExpansionSharedPtr
NekDouble Nektar::MultiRegions::ExpList3DHomogeneous1D::v_L2 ( const Array< OneD, const NekDouble > &  inarray,
const Array< OneD, const NekDouble > &  soln = NullNekDouble1DArray 
)
protectedvirtual

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

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

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

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 478 of file ExpList3DHomogeneous1D.cpp.

References Nektar::MultiRegions::ExpList::m_comm, Nektar::MultiRegions::ExpListHomogeneous1D::m_homogeneousBasis, Nektar::MultiRegions::ExpListHomogeneous1D::m_lhom, Nektar::MultiRegions::ExpListHomogeneous1D::m_planes, Nektar::MultiRegions::ExpListHomogeneous1D::m_transposition, Nektar::NullNekDouble1DArray, and Nektar::LibUtilities::ReduceSum.

481  {
482  int cnt = 0;
483  NekDouble errL2, err = 0.0;
484  Array<OneD, const NekDouble> w = m_homogeneousBasis->GetW();
485  Array<OneD, NekDouble> local_w(m_planes.num_elements());
486 
487  for(int n = 0; n < m_planes.num_elements(); n++)
488  {
489  local_w[n] = w[m_transposition->GetPlaneID(n)];
490  }
491 
492  if (soln == NullNekDouble1DArray)
493  {
494  for(int n = 0; n < m_planes.num_elements(); ++n)
495  {
496  errL2 = m_planes[n]->L2(inarray + cnt);
497  cnt += m_planes[n]->GetTotPoints();
498  err += errL2*errL2*local_w[n]*m_lhom*0.5;
499  }
500  }
501  else
502  {
503  for(int n = 0; n < m_planes.num_elements(); ++n)
504  {
505  errL2 = m_planes[n]->L2(inarray + cnt, soln + cnt);
506  cnt += m_planes[n]->GetTotPoints();
507  err += errL2*errL2*local_w[n]*m_lhom*0.5;
508  }
509  }
510  m_comm->GetColumnComm()->AllReduce(err, LibUtilities::ReduceSum);
511 
512  return sqrt(err);
513  }
static Array< OneD, NekDouble > NullNekDouble1DArray
LibUtilities::TranspositionSharedPtr m_transposition
NekDouble m_lhom
Width of homogeneous direction.
LibUtilities::BasisSharedPtr m_homogeneousBasis
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
Array< OneD, ExpListSharedPtr > m_planes
double NekDouble
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: ExpList.h:907
void Nektar::MultiRegions::ExpList3DHomogeneous1D::v_WriteTecplotConnectivity ( std::ostream &  outfile,
int  expansion 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 321 of file ExpList3DHomogeneous1D.cpp.

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

322  {
323  ASSERTL0(expansion == -1, "Multi-zone output not supported for homogeneous expansions.");
324 
325  const int nPtsPlane = m_planes[0]->GetNpoints();
326  const int nElmt = m_planes[0]->GetExpSize();
327  const int nPlanes = m_planes.num_elements();
328 
329  int cnt = 0;
330  int cnt2 = 0;
331  for (int i = 0; i < nElmt; ++i)
332  {
333  const int np0 = (*m_exp)[i]->GetNumPoints(0);
334  const int np1 = (*m_exp)[i]->GetNumPoints(1);
335 
336  for (int n = 1; n < nPlanes; ++n)
337  {
338  const int o1 = (n-1) * nPtsPlane;
339  const int o2 = n * nPtsPlane;
340  for (int j = 1; j < np1; ++j)
341  {
342  for(int k = 1; k < np0; ++k)
343  {
344  outfile << cnt + (j-1)*np0 + (k-1) + o1 + 1 << " ";
345  outfile << cnt + (j-1)*np0 + (k-1) + o2 + 1 << " ";
346  outfile << cnt + (j-1)*np0 + k + o2 + 1 << " ";
347  outfile << cnt + (j-1)*np0 + k + o1 + 1 << " ";
348  outfile << cnt + j *np0 + (k-1) + o1 + 1 << " ";
349  outfile << cnt + j *np0 + (k-1) + o2 + 1 << " ";
350  outfile << cnt + j *np0 + k + o2 + 1 << " ";
351  outfile << cnt + j *np0 + k + o1 + 1 << endl;
352  cnt2++;
353  }
354  }
355  }
356 
357  cnt += np0*np1;
358  }
359  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
Array< OneD, ExpListSharedPtr > m_planes
void Nektar::MultiRegions::ExpList3DHomogeneous1D::v_WriteVtkPieceHeader ( std::ostream &  outfile,
int  expansion,
int  istrip 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 362 of file ExpList3DHomogeneous1D.cpp.

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

366  {
367  // If there is only one plane (e.g. HalfMode), we write a 2D plane.
368  if (m_planes.num_elements() == 1)
369  {
370  m_planes[0]->WriteVtkPieceHeader(outfile, expansion);
371  return;
372  }
373 
374  // If we are using Fourier points, output extra plane to fill domain
375  int outputExtraPlane = 0;
376  if ( m_homogeneousBasis->GetBasisType() == LibUtilities::eFourier
377  && m_homogeneousBasis->GetPointsType() ==
379  {
380  outputExtraPlane = 1;
381  }
382 
383  int i,j,k;
384  int nq0 = (*m_exp)[expansion]->GetNumPoints(0);
385  int nq1 = (*m_exp)[expansion]->GetNumPoints(1);
386  int nq2 = m_planes.num_elements() + outputExtraPlane;
387  int ntot = nq0*nq1*nq2;
388  int ntotminus = (nq0-1)*(nq1-1)*(nq2-1);
389 
390  Array<OneD,NekDouble> coords[3];
391  coords[0] = Array<OneD,NekDouble>(ntot);
392  coords[1] = Array<OneD,NekDouble>(ntot);
393  coords[2] = Array<OneD,NekDouble>(ntot);
394  GetCoords(expansion,coords[0],coords[1],coords[2]);
395 
396  if (outputExtraPlane)
397  {
398  // Copy coords[0] and coords[1] to extra plane
399  Array<OneD,NekDouble> tmp;
400  Vmath::Vcopy (nq0*nq1, coords[0], 1,
401  tmp = coords[0] + (nq2-1)*nq0*nq1, 1);
402  Vmath::Vcopy (nq0*nq1, coords[1], 1,
403  tmp = coords[1] + (nq2-1)*nq0*nq1, 1);
404  // Fill coords[2] for extra plane
405  NekDouble z = coords[2][nq0*nq1*m_planes.num_elements()-1] +
406  (coords[2][nq0*nq1] - coords[2][0]);
407  Vmath::Fill(nq0*nq1, z, tmp = coords[2] + (nq2-1)*nq0*nq1, 1);
408  }
409 
410  NekDouble DistStrip;
411  m_session->LoadParameter("DistStrip", DistStrip, 0);
412  // Reset the z-coords for homostrips
413  for(int i = 0; i < ntot; i++)
414  {
415  coords[2][i] += istrip*DistStrip;
416  }
417 
418  outfile << " <Piece NumberOfPoints=\""
419  << ntot << "\" NumberOfCells=\""
420  << ntotminus << "\">" << endl;
421  outfile << " <Points>" << endl;
422  outfile << " <DataArray type=\"Float64\" "
423  << "NumberOfComponents=\"3\" format=\"ascii\">" << endl;
424  outfile << " ";
425  for (i = 0; i < ntot; ++i)
426  {
427  for (j = 0; j < 3; ++j)
428  {
429  outfile << coords[j][i] << " ";
430  }
431  outfile << endl;
432  }
433  outfile << endl;
434  outfile << " </DataArray>" << endl;
435  outfile << " </Points>" << endl;
436  outfile << " <Cells>" << endl;
437  outfile << " <DataArray type=\"Int32\" "
438  << "Name=\"connectivity\" format=\"ascii\">" << endl;
439  for (i = 0; i < nq0-1; ++i)
440  {
441  for (j = 0; j < nq1-1; ++j)
442  {
443  for (k = 0; k < nq2-1; ++k)
444  {
445  outfile << k*nq0*nq1 + j*nq0 + i << " "
446  << k*nq0*nq1 + j*nq0 + i + 1 << " "
447  << k*nq0*nq1 + (j+1)*nq0 + i + 1 << " "
448  << k*nq0*nq1 + (j+1)*nq0 + i << " "
449  << (k+1)*nq0*nq1 + j*nq0 + i << " "
450  << (k+1)*nq0*nq1 + j*nq0 + i + 1 << " "
451  << (k+1)*nq0*nq1 + (j+1)*nq0 + i + 1 << " "
452  << (k+1)*nq0*nq1 + (j+1)*nq0 + i << endl;
453  }
454  }
455  }
456  outfile << endl;
457  outfile << " </DataArray>" << endl;
458  outfile << " <DataArray type=\"Int32\" "
459  << "Name=\"offsets\" format=\"ascii\">" << endl;
460  for (i = 0; i < ntotminus; ++i)
461  {
462  outfile << i*8+8 << " ";
463  }
464  outfile << endl;
465  outfile << " </DataArray>" << endl;
466  outfile << " <DataArray type=\"UInt8\" "
467  << "Name=\"types\" format=\"ascii\">" << endl;
468  for (i = 0; i < ntotminus; ++i)
469  {
470  outfile << "12 ";
471  }
472  outfile << endl;
473  outfile << " </DataArray>" << endl;
474  outfile << " </Cells>" << endl;
475  outfile << " <PointData>" << endl;
476  }
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:46
Fourier Expansion .
Definition: BasisType.h:52
LibUtilities::BasisSharedPtr m_homogeneousBasis
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:64
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:910
Array< OneD, ExpListSharedPtr > m_planes
double NekDouble
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 . ...
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047