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

#include <DisContField3DHomogeneous2D.h>

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

Public Member Functions

 DisContField3DHomogeneous2D ()
 DisContField3DHomogeneous2D (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::BasisKey &HomoBasis_y, const LibUtilities::BasisKey &HomoBasis_z, const NekDouble lhom_y, const NekDouble lhom_z, const bool useFFT, const bool dealiasing)
 DisContField3DHomogeneous2D (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::BasisKey &HomoBasis_y, const LibUtilities::BasisKey &HomoBasis_z, const NekDouble lhom_y, const NekDouble lhom_z, const bool useFFT, const bool dealiasing, const SpatialDomains::MeshGraphSharedPtr &graph1D, const std::string &variable)
 DisContField3DHomogeneous2D (const DisContField3DHomogeneous2D &In, const bool DeclareLinesSetCoeffPhys=true)
 Copy constructor.
virtual ~DisContField3DHomogeneous2D ()
 Destructor.
void SetupBoundaryConditions (const LibUtilities::BasisKey &HomoBasis_y, const LibUtilities::BasisKey &HomoBasis_z, const NekDouble lhom_y, const NekDouble lhom_z, SpatialDomains::BoundaryConditions &bcs)
void EvaluateBoundaryConditions (const NekDouble time=0.0, const std::string varName="")
const Array< OneD, const
MultiRegions::ExpListSharedPtr > & 
GetBndCondExpansions ()
const Array< OneD, const
SpatialDomains::BoundaryConditionShPtr > & 
GetBndConditions ()
boost::shared_ptr< ExpList > & UpdateBndCondExpansion (int i)
Array< OneD,
SpatialDomains::BoundaryConditionShPtr > & 
UpdateBndConditions ()
void GetBoundaryToElmtMap (Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
 Set up a list of element ids and edge ids the link to the boundary conditions.
- Public Member Functions inherited from Nektar::MultiRegions::ExpList3DHomogeneous2D
 ExpList3DHomogeneous2D ()
 Default constructor.
 ExpList3DHomogeneous2D (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::BasisKey &HomoBasis_y, const LibUtilities::BasisKey &HomoBasis_z, const NekDouble lhom_y, const NekDouble lhom_z, const bool useFFT, const bool dealiasing)
 ExpList3DHomogeneous2D (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::BasisKey &HomoBasis_y, const LibUtilities::BasisKey &HomoBasis_z, const NekDouble lhom_y, const NekDouble lhom_z, const bool useFFT, const bool dealiasing, const SpatialDomains::MeshGraphSharedPtr &graph1D)
 Sets up a list of local expansions based on an input mesh.
 ExpList3DHomogeneous2D (const ExpList3DHomogeneous2D &In, const bool DeclareLinesSetCoeffPhys=true)
 Copy constructor.
virtual ~ExpList3DHomogeneous2D ()
 Destructor.
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 GetCoords (const int eid, Array< OneD, NekDouble > &xc0, Array< OneD, NekDouble > &xc1, Array< OneD, NekDouble > &xc2)
- Public Member Functions inherited from Nektar::MultiRegions::ExpListHomogeneous2D
 ExpListHomogeneous2D ()
 Default constructor.
 ExpListHomogeneous2D (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::BasisKey &HomoBasis_y, const LibUtilities::BasisKey &HomoBasis_z, const NekDouble ly, const NekDouble lz, const bool useFFT, const bool dealiasing)
 ExpListHomogeneous2D (const ExpListHomogeneous2D &In)
 Copy constructor.
virtual ~ExpListHomogeneous2D ()
 Destructor.
void Homogeneous2DTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool IsForwards, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
void HomogeneousFwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
void HomogeneousBwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
void DealiasedProd (const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
void SetPaddingBase (void)
void PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
 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)
- 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 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)
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)
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 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 ()
LibUtilities::BasisSharedPtr GetHomogeneousBasis (void)
boost::shared_ptr< ExpList > & GetPlane (int n)

Public Attributes

Array< OneD, int > m_BCtoElmMap
 Storage space for the boundary to element and boundary to trace map. This member variable is really allocated just in case a boundary expansion recasting is required at the solver level. Otherwise is the 2 vectors are not filled up. If is needed all the funcitons whihc require to use this map do not have to recalculate it anymore.
Array< OneD, int > m_BCtoEdgMap

Protected Member Functions

virtual void v_GetBoundaryToElmtMap (Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
virtual map< int,
RobinBCInfoSharedPtr
v_GetRobinBCInfo ()
- Protected Member Functions inherited from Nektar::MultiRegions::ExpList3DHomogeneous2D
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.
virtual void v_GetCoords (Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1, Array< OneD, NekDouble > &coord_2)
virtual void v_WriteTecplotZone (std::ofstream &outfile, int expansion)
virtual void v_WriteVtkPieceHeader (std::ofstream &outfile, int expansion)
virtual NekDouble v_L2 (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
- Protected Member Functions inherited from Nektar::MultiRegions::ExpListHomogeneous2D
DNekBlkMatSharedPtr GenHomogeneous2DBlockMatrix (Homogeneous2DMatType mattype, CoeffState coeffstate=eLocal) const
DNekBlkMatSharedPtr GetHomogeneous2DBlockMatrix (Homogeneous2DMatType mattype, CoeffState coeffstate=eLocal) const
virtual int v_GetNumElmts (void)
virtual void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
virtual void v_FwdTrans_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
virtual void v_BwdTrans_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
virtual void v_IProductWRTBase_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual std::vector
< LibUtilities::FieldDefinitionsSharedPtr
v_GetFieldDefinitions (void)
virtual void v_GetFieldDefinitions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef)
virtual void v_AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
virtual void v_AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, Array< OneD, NekDouble > &coeffs)
virtual void v_ExtractDataToCoeffs (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
 Extract data from raw field data into expansion list.
virtual void v_HomogeneousFwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
virtual void v_HomogeneousBwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
virtual void v_DealiasedProd (const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
virtual void v_PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
virtual void v_PhysDeriv (Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
- Protected Member Functions inherited from Nektar::MultiRegions::ExpList
boost::shared_ptr< DNekMatGenGlobalMatrixFull (const GlobalLinSysKey &mkey, const boost::shared_ptr< AssemblyMapCG > &locToGloMap)
const DNekScalBlkMatSharedPtr GenBlockMatrix (const GlobalMatrixKey &gkey)
 This function assembles the block diagonal matrix of local matrices of the type mtype.
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_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_ReadGlobalOptimizationParameters ()
virtual void v_ExtractCoeffsToCoeffs (const boost::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
virtual void v_WriteTecplotHeader (std::ofstream &outfile, std::string var="")
virtual void v_WriteTecplotField (std::ofstream &outfile, int expansion)
virtual void v_WriteTecplotConnectivity (std::ofstream &outfile, int expansion)
virtual void v_WriteVtkPieceData (std::ofstream &outfile, int expansion, std::string var)
virtual NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray)
virtual Array< OneD, const
NekDouble
v_HomogeneousEnergy (void)
virtual
LibUtilities::TranspositionSharedPtr 
v_GetTransposition (void)
virtual NekDouble v_GetHomoLen (void)
virtual Array< OneD, const
unsigned int > 
v_GetZIDs (void)
virtual Array< OneD, const
unsigned int > 
v_GetYIDs (void)
virtual void v_PhysInterp1DScaled (const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_PhysGalerkinProjection1DScaled (const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void ExtractFileBCs (const std::string &fileName, const std::string &varName, const boost::shared_ptr< ExpList > locExpList)

Protected Attributes

Array< OneD,
MultiRegions::ExpListSharedPtr
m_bndCondExpansions
Array< OneD,
SpatialDomains::BoundaryConditionShPtr
m_bndConditions

Private Member Functions

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_EvaluateBoundaryConditions (const NekDouble time=0.0, const std::string varName="", const NekDouble x2_in=NekConstants::kNekUnsetDouble, const NekDouble x3_in=NekConstants::kNekUnsetDouble)
virtual const Array< OneD,
const boost::shared_ptr
< ExpList > > & 
v_GetBndCondExpansions (void)
virtual const Array< OneD,
const
SpatialDomains::BoundaryConditionShPtr > & 
v_GetBndConditions ()
virtual boost::shared_ptr
< ExpList > & 
v_UpdateBndCondExpansion (int i)
virtual Array< OneD,
SpatialDomains::BoundaryConditionShPtr > & 
v_UpdateBndConditions ()

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

Definition at line 47 of file DisContField3DHomogeneous2D.h.

Constructor & Destructor Documentation

Nektar::MultiRegions::DisContField3DHomogeneous2D::DisContField3DHomogeneous2D ( void  )
Nektar::MultiRegions::DisContField3DHomogeneous2D::DisContField3DHomogeneous2D ( const LibUtilities::SessionReaderSharedPtr pSession,
const LibUtilities::BasisKey HomoBasis_y,
const LibUtilities::BasisKey HomoBasis_z,
const NekDouble  lhom_y,
const NekDouble  lhom_z,
const bool  useFFT,
const bool  dealiasing 
)

Definition at line 53 of file DisContField3DHomogeneous2D.cpp.

:
ExpList3DHomogeneous2D(pSession,HomoBasis_y,HomoBasis_z,lhom_y,lhom_z,useFFT,dealiasing),
{
}
Nektar::MultiRegions::DisContField3DHomogeneous2D::DisContField3DHomogeneous2D ( const LibUtilities::SessionReaderSharedPtr pSession,
const LibUtilities::BasisKey HomoBasis_y,
const LibUtilities::BasisKey HomoBasis_z,
const NekDouble  lhom_y,
const NekDouble  lhom_z,
const bool  useFFT,
const bool  dealiasing,
const SpatialDomains::MeshGraphSharedPtr graph1D,
const std::string &  variable 
)

Definition at line 84 of file DisContField3DHomogeneous2D.cpp.

References Nektar::MultiRegions::ExpList::GetExp(), Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::MultiRegions::ExpList::m_exp, Nektar::MultiRegions::ExpList::m_globalOptParam, Nektar::MultiRegions::ExpListHomogeneous2D::m_homogeneousBasis_y, Nektar::MultiRegions::ExpListHomogeneous2D::m_homogeneousBasis_z, Nektar::MultiRegions::ExpListHomogeneous2D::m_lines, Nektar::MultiRegions::ExpList3DHomogeneous2D::SetCoeffPhys(), and SetupBoundaryConditions().

:
ExpList3DHomogeneous2D(pSession,HomoBasis_y,HomoBasis_z,lhom_y,lhom_z,useFFT,dealiasing),
{
int i,n,nel;
SpatialDomains::BoundaryConditions bcs(pSession, graph1D);
//
m_lines[0] = line_zero = MemoryManager<DisContField1D>::AllocateSharedPtr(pSession,graph1D,variable);
nel = m_lines[0]->GetExpSize();
for(i = 0; i < nel; ++i)
{
(*m_exp).push_back(m_lines[0]->GetExp(i));
}
int nylines = m_homogeneousBasis_y->GetNumPoints();
int nzlines = m_homogeneousBasis_z->GetNumPoints();
for(n = 1; n < nylines*nzlines; ++n)
{
for(i = 0; i < nel; ++i)
{
(*m_exp).push_back((*m_exp)[i]);
}
}
// Setup Default optimisation information.
nel = GetExpSize();
SetupBoundaryConditions(HomoBasis_y,HomoBasis_z,lhom_y,lhom_z,bcs);
}
Nektar::MultiRegions::DisContField3DHomogeneous2D::DisContField3DHomogeneous2D ( const DisContField3DHomogeneous2D In,
const bool  DeclareLinesSetCoeffPhys = true 
)

Copy constructor.

Definition at line 66 of file DisContField3DHomogeneous2D.cpp.

References Nektar::MultiRegions::ExpListHomogeneous2D::m_lines, and Nektar::MultiRegions::ExpList3DHomogeneous2D::SetCoeffPhys().

:
m_bndCondExpansions (In.m_bndCondExpansions),
m_bndConditions (In.m_bndConditions)
{
if(DeclareLinesSetCoeffPhys)
{
DisContField1DSharedPtr zero_line = boost::dynamic_pointer_cast<DisContField1D> (In.m_lines[0]);
for(int n = 0; n < m_lines.num_elements(); ++n)
{
}
}
}
Nektar::MultiRegions::DisContField3DHomogeneous2D::~DisContField3DHomogeneous2D ( )
virtual

Destructor.

Definition at line 135 of file DisContField3DHomogeneous2D.cpp.

{
}

Member Function Documentation

void Nektar::MultiRegions::DisContField3DHomogeneous2D::EvaluateBoundaryConditions ( const NekDouble  time = 0.0,
const std::string  varName = "" 
)

Definition at line 178 of file DisContField3DHomogeneous2D.cpp.

References Nektar::SpatialDomains::eTimeDependent, Nektar::MultiRegions::ExpList::GetCoeffs(), m_bndCondExpansions, m_bndConditions, Nektar::MultiRegions::ExpListHomogeneous2D::m_homogeneousBasis_y, Nektar::MultiRegions::ExpListHomogeneous2D::m_homogeneousBasis_z, Nektar::MultiRegions::ExpListHomogeneous2D::m_lhom_y, Nektar::MultiRegions::ExpListHomogeneous2D::m_lhom_z, Nektar::MultiRegions::ExpListHomogeneous2D::m_lines, Nektar::MultiRegions::ExpListHomogeneous2D::m_ny, Nektar::MultiRegions::ExpListHomogeneous2D::m_nz, and Nektar::MultiRegions::ExpList::UpdateCoeffs().

Referenced by SetupBoundaryConditions(), and v_EvaluateBoundaryConditions().

{
int n, m;
const Array<OneD, const NekDouble> y = m_homogeneousBasis_y->GetZ();
const Array<OneD, const NekDouble> z = m_homogeneousBasis_z->GetZ();
for (n = 0; n < m_nz; ++n)
{
for (m = 0; m < m_ny; ++m)
{
time, varName, 0.5*m_lhom_y*(1.0+y[m]),
0.5*m_lhom_z*(1.0+z[n]));
}
}
// Fourier transform coefficient space boundary values
for (n = 0; n < m_bndCondExpansions.num_elements(); ++n)
{
if (time == 0.0 || m_bndConditions[n]->GetUserDefined() ==
{
m_bndCondExpansions[n]->HomogeneousFwdTrans(
}
}
}
const Array< OneD, const MultiRegions::ExpListSharedPtr > & Nektar::MultiRegions::DisContField3DHomogeneous2D::GetBndCondExpansions ( )
inline

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 153 of file DisContField3DHomogeneous2D.h.

References m_bndCondExpansions.

Referenced by v_GetBndCondExpansions().

{
}
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & Nektar::MultiRegions::DisContField3DHomogeneous2D::GetBndConditions ( )
inline

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 158 of file DisContField3DHomogeneous2D.h.

References m_bndConditions.

Referenced by v_GetBndConditions().

{
}
void Nektar::MultiRegions::DisContField3DHomogeneous2D::GetBoundaryToElmtMap ( Array< OneD, int > &  ElmtID,
Array< OneD, int > &  EdgeID 
)

Set up a list of element ids and edge ids the link to the boundary conditions.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 287 of file DisContField3DHomogeneous2D.cpp.

References m_BCtoEdgMap, m_BCtoElmMap, Nektar::MultiRegions::ExpListHomogeneous2D::m_lines, and Vmath::Vcopy().

Referenced by v_GetBoundaryToElmtMap().

{
if(m_BCtoElmMap.num_elements() == 0)
{
Array<OneD, int> ElmtID_tmp;
Array<OneD, int> EdgeID_tmp;
m_lines[0]->GetBoundaryToElmtMap(ElmtID_tmp,EdgeID_tmp);
int nel_per_lines = m_lines[0]->GetExpSize();
int nlines = m_lines.num_elements();
int MapSize = ElmtID_tmp.num_elements();
ElmtID = Array<OneD, int>(nlines*MapSize);
EdgeID = Array<OneD, int>(nlines*MapSize);
for(int i = 0; i < nlines; i++)
{
for(int j = 0; j < MapSize; j++)
{
ElmtID[j+i*MapSize] = ElmtID_tmp[j] + i*nel_per_lines;
EdgeID[j+i*MapSize] = EdgeID_tmp[j];
}
}
m_BCtoElmMap = Array<OneD, int>(nlines*MapSize);
m_BCtoEdgMap = Array<OneD, int>(nlines*MapSize);
Vmath::Vcopy(nlines*MapSize,ElmtID,1,m_BCtoElmMap,1);
Vmath::Vcopy(nlines*MapSize,EdgeID,1,m_BCtoEdgMap,1);
}
else
{
int MapSize = m_BCtoElmMap.num_elements();
ElmtID = Array<OneD, int>(MapSize);
EdgeID = Array<OneD, int>(MapSize);
Vmath::Vcopy(MapSize,m_BCtoElmMap,1,ElmtID,1);
Vmath::Vcopy(MapSize,m_BCtoEdgMap,1,EdgeID,1);
}
}
void Nektar::MultiRegions::DisContField3DHomogeneous2D::SetupBoundaryConditions ( const LibUtilities::BasisKey HomoBasis_y,
const LibUtilities::BasisKey HomoBasis_z,
const NekDouble  lhom_y,
const NekDouble  lhom_z,
SpatialDomains::BoundaryConditions bcs 
)

Definition at line 140 of file DisContField3DHomogeneous2D.cpp.

References EvaluateBoundaryConditions(), Nektar::SpatialDomains::BoundaryConditions::GetBoundaryRegions(), m_bndCondExpansions, m_bndConditions, Nektar::MultiRegions::ExpListHomogeneous2D::m_lines, Nektar::MultiRegions::ExpList::m_session, and Nektar::MultiRegions::ExpListHomogeneous2D::m_useFFT.

Referenced by Nektar::MultiRegions::ContField3DHomogeneous2D::ContField3DHomogeneous2D(), and DisContField3DHomogeneous2D().

{
int i,n;
// Setup an ExpList1DHomogeneous2D expansion for boundary
// conditions and link to class declared in m_lines.
int nlines = m_lines.num_elements();
int nbnd = bregions.size();
m_bndCondExpansions = Array<OneD,MultiRegions::ExpListSharedPtr>(nbnd);
Array<OneD, MultiRegions::ExpListSharedPtr> LinesBndCondExp(nlines);
m_bndConditions = m_lines[0]->UpdateBndConditions();
for(i = 0; i < nbnd; ++i)
{
for(n = 0; n < nlines; ++n)
{
LinesBndCondExp[n] = m_lines[n]->UpdateBndCondExpansion(i);
}
m_bndCondExpansions[i] = MemoryManager<ExpList1DHomogeneous2D>::AllocateSharedPtr(m_session,HomoBasis_y,HomoBasis_z,lhom_y,lhom_z,m_useFFT,false,LinesBndCondExp);
}
}
MultiRegions::ExpListSharedPtr & Nektar::MultiRegions::DisContField3DHomogeneous2D::UpdateBndCondExpansion ( int  i)
inline

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 163 of file DisContField3DHomogeneous2D.h.

References m_bndCondExpansions.

Referenced by v_UpdateBndCondExpansion().

{
}
Array< OneD, SpatialDomains::BoundaryConditionShPtr > & Nektar::MultiRegions::DisContField3DHomogeneous2D::UpdateBndConditions ( )
inline

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 168 of file DisContField3DHomogeneous2D.h.

References m_bndConditions.

Referenced by v_UpdateBndConditions().

{
}
void Nektar::MultiRegions::DisContField3DHomogeneous2D::v_EvaluateBoundaryConditions ( const NekDouble  time = 0.0,
const std::string  varName = "",
const NekDouble  x2_in = NekConstants::kNekUnsetDouble,
const NekDouble  x3_in = NekConstants::kNekUnsetDouble 
)
privatevirtual

Definition at line 258 of file DisContField3DHomogeneous2D.cpp.

References EvaluateBoundaryConditions().

{
}
const Array< OneD, const boost::shared_ptr< ExpList > > & Nektar::MultiRegions::DisContField3DHomogeneous2D::v_GetBndCondExpansions ( void  )
privatevirtual

Definition at line 267 of file DisContField3DHomogeneous2D.cpp.

References GetBndCondExpansions().

{
}
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & Nektar::MultiRegions::DisContField3DHomogeneous2D::v_GetBndConditions ( void  )
privatevirtual

Definition at line 272 of file DisContField3DHomogeneous2D.cpp.

References GetBndConditions().

{
return GetBndConditions();
}
virtual void Nektar::MultiRegions::DisContField3DHomogeneous2D::v_GetBoundaryToElmtMap ( Array< OneD, int > &  ElmtID,
Array< OneD, int > &  EdgeID 
)
inlineprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 113 of file DisContField3DHomogeneous2D.h.

References GetBoundaryToElmtMap().

{
GetBoundaryToElmtMap(ElmtID,EdgeID);
}
virtual map<int, RobinBCInfoSharedPtr> Nektar::MultiRegions::DisContField3DHomogeneous2D::v_GetRobinBCInfo ( void  )
inlineprotectedvirtual
Todo:
Fix Robin BCs for homogeneous case

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 120 of file DisContField3DHomogeneous2D.h.

{
return map<int, RobinBCInfoSharedPtr>();
}
void Nektar::MultiRegions::DisContField3DHomogeneous2D::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 
)
privatevirtual

Definition at line 209 of file DisContField3DHomogeneous2D.cpp.

References Nektar::StdRegions::eFactorLambda, Nektar::MultiRegions::ExpListHomogeneous2D::HomogeneousFwdTrans(), Nektar::MultiRegions::ExpListHomogeneous2D::m_homogeneousBasis_y, Nektar::MultiRegions::ExpListHomogeneous2D::m_homogeneousBasis_z, Nektar::MultiRegions::ExpListHomogeneous2D::m_lhom_y, Nektar::MultiRegions::ExpListHomogeneous2D::m_lhom_z, Nektar::MultiRegions::ExpListHomogeneous2D::m_lines, and Nektar::MultiRegions::ExpList::m_WaveSpace.

{
int n,m;
int cnt = 0;
int cnt1 = 0;
int nhom_modes_y = m_homogeneousBasis_y->GetNumModes();
int nhom_modes_z = m_homogeneousBasis_z->GetNumModes();
NekDouble beta_y;
NekDouble beta_z;
Array<OneD, NekDouble> e_out;
Array<OneD, NekDouble> fce(inarray.num_elements());
// Fourier transform forcing function
{
fce = inarray;
}
else
{
HomogeneousFwdTrans(inarray,fce);
}
for(n = 0; n < nhom_modes_z; ++n)
{
for(m = 0; m < nhom_modes_y; ++m)
{
beta_z = 2*M_PI*(n/2)/m_lhom_z;
beta_y = 2*M_PI*(m/2)/m_lhom_y;
new_factors = factors;
new_factors[StdRegions::eFactorLambda] += beta_y*beta_y + beta_z*beta_z;
m_lines[n]->HelmSolve(fce + cnt,
e_out = outarray + cnt1,
flags, new_factors, varcoeff, dirForcing);
cnt += m_lines[n]->GetTotPoints();
cnt1 += m_lines[n]->GetNcoeffs();
}
}
}
boost::shared_ptr< ExpList > & Nektar::MultiRegions::DisContField3DHomogeneous2D::v_UpdateBndCondExpansion ( int  i)
privatevirtual

Definition at line 277 of file DisContField3DHomogeneous2D.cpp.

References UpdateBndCondExpansion().

{
}
Array< OneD, SpatialDomains::BoundaryConditionShPtr > & Nektar::MultiRegions::DisContField3DHomogeneous2D::v_UpdateBndConditions ( )
privatevirtual

Definition at line 282 of file DisContField3DHomogeneous2D.cpp.

References UpdateBndConditions().

{
}

Member Data Documentation

Array<OneD, int> Nektar::MultiRegions::DisContField3DHomogeneous2D::m_BCtoEdgMap

Definition at line 105 of file DisContField3DHomogeneous2D.h.

Referenced by GetBoundaryToElmtMap().

Array<OneD, int> Nektar::MultiRegions::DisContField3DHomogeneous2D::m_BCtoElmMap

Storage space for the boundary to element and boundary to trace map. This member variable is really allocated just in case a boundary expansion recasting is required at the solver level. Otherwise is the 2 vectors are not filled up. If is needed all the funcitons whihc require to use this map do not have to recalculate it anymore.

Definition at line 104 of file DisContField3DHomogeneous2D.h.

Referenced by GetBoundaryToElmtMap().

Array<OneD,MultiRegions::ExpListSharedPtr> Nektar::MultiRegions::DisContField3DHomogeneous2D::m_bndCondExpansions
protected

Definition at line 109 of file DisContField3DHomogeneous2D.h.

Referenced by EvaluateBoundaryConditions(), GetBndCondExpansions(), SetupBoundaryConditions(), and UpdateBndCondExpansion().

Array<OneD,SpatialDomains::BoundaryConditionShPtr> Nektar::MultiRegions::DisContField3DHomogeneous2D::m_bndConditions
protected

Definition at line 111 of file DisContField3DHomogeneous2D.h.

Referenced by EvaluateBoundaryConditions(), GetBndConditions(), SetupBoundaryConditions(), and UpdateBndConditions().