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

#include <ContField3D.h>

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

Public Member Functions

 ContField3D ()
 ContField3D (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &graph3D, const std::string &variable="DefaultVar", const bool CheckIfSingularSystem=false)
 Construct a global continuous field.
 ContField3D (const ContField3D &In, const SpatialDomains::MeshGraphSharedPtr &graph3D, const std::string &variable, const bool CheckIfSingularSystem=false)
 Construct a global continuous field with solution type based on another field but using a separate input mesh and boundary conditions.
 ContField3D (const ContField3D &In)
virtual ~ContField3D ()
const Array< OneD, const
MultiRegions::ExpListSharedPtr > & 
GetBndCondExpansions ()
const Array< OneD, const
MultiRegions::ExpListSharedPtr > & 
GetBndCondExp ()
 This function return the boundary conditions expansion.
void GenerateDirBndCondForcing (const GlobalLinSysKey &key, Array< OneD, NekDouble > &inout, Array< OneD, NekDouble > &outarray)
void GlobalToLocal (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void LocalToGlobal (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void Assemble ()
void Assemble (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
const AssemblyMapCGSharedPtrGetLocalToGlobalMap () const
int GetGlobalMatrixNnz (const GlobalMatrixKey &gkey)
- Public Member Functions inherited from Nektar::MultiRegions::DisContField3D
 DisContField3D ()
 Default constructor.
 DisContField3D (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &graph3D, const std::string &variable, const bool SetUpJustDG=true)
 Constructs a global discontinuous field based on an input mesh with boundary conditions.
 DisContField3D (const DisContField3D &In, const SpatialDomains::MeshGraphSharedPtr &graph3D, const std::string &variable, const bool SetUpJustDG=false)
 DisContField3D (const DisContField3D &In)
 Constructs a global discontinuous field based on another discontinuous field.
virtual ~DisContField3D ()
 Destructor.
GlobalLinSysSharedPtr GetGlobalBndLinSys (const GlobalLinSysKey &mkey)
void EvaluateHDGPostProcessing (Array< OneD, NekDouble > &outarray)
 Evaluate HDG post-processing to increase polynomial order of solution.
- Public Member Functions inherited from Nektar::MultiRegions::ExpList3D
 ExpList3D ()
 Default constructor.
 ExpList3D (const ExpList3D &In)
 Copy constructor.
 ExpList3D (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::BasisKey &TBa, const LibUtilities::BasisKey &TBb, const LibUtilities::BasisKey &TBc, const LibUtilities::BasisKey &HBa, const LibUtilities::BasisKey &HBb, const LibUtilities::BasisKey &HBc, const SpatialDomains::MeshGraphSharedPtr &graph3D, const LibUtilities::PointsType TetNb=LibUtilities::SIZE_PointsType)
 ExpList3D (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &graph3D, const std::string &variable="DefaultVar")
 Sets up a list of local expansions based on an input mesh.
 ExpList3D (const SpatialDomains::ExpansionMap &expansions)
 Sets up a list of local expansions based on an expansion vector.
virtual ~ExpList3D ()
 Destructor.
- Public Member Functions inherited from Nektar::MultiRegions::ExpList
 ExpList ()
 The default constructor.
 ExpList (const LibUtilities::SessionReaderSharedPtr &pSession)
 The default constructor.
 ExpList (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 The default constructor.
 ExpList (const ExpList &in, const bool DeclareCoeffPhysArrays=true)
 The copy constructor.
virtual ~ExpList ()
 The default destructor.
int GetNcoeffs (void) const
 Returns the total number of local degrees of freedom $N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_m$.
int GetNcoeffs (const int eid) const
 Returns the total number of local degrees of freedom for element eid.
ExpansionType GetExpType (void)
 Returns the type of the expansion.
void SetExpType (ExpansionType Type)
 Returns the type of the expansion.
int EvalBasisNumModesMax (void) const
 Evaulates the maximum number of modes in the elemental basis order over all elements.
const Array< OneD, int > EvalBasisNumModesMaxPerExp (void) const
 Returns the vector of the number of modes in the elemental basis order over all elements.
int GetTotPoints (void) const
 Returns the total number of quadrature points m_npoints $=Q_{\mathrm{tot}}$.
int GetTotPoints (const int eid) const
 Returns the total number of quadrature points for eid's element $=Q_{\mathrm{tot}}$.
int GetNpoints (void) const
 Returns the total number of quadrature points m_npoints $=Q_{\mathrm{tot}}$.
int Get1DScaledTotPoints (const NekDouble scale) const
 Returns the total number of qudature points scaled by the factor scale on each 1D direction.
void SetWaveSpace (const bool wavespace)
 Sets the wave space to the one of the possible configuration true or false.
void SetModifiedBasis (const bool modbasis)
 Set Modified Basis for the stability analysis.
void SetPhys (int i, NekDouble val)
 Set the i th value of m_phys to value val.
bool GetWaveSpace (void) const
 This function returns the third direction expansion condition, which can be in wave space (coefficient) or not It is stored in the variable m_WaveSpace.
void SetPhys (const Array< OneD, const NekDouble > &inarray)
 Fills the array m_phys.
void SetPhysArray (Array< OneD, NekDouble > &inarray)
 Sets the array m_phys.
void SetPhysState (const bool physState)
 This function manually sets whether the array of physical values $\boldsymbol{u}_l$ (implemented as m_phys) is filled or not.
bool GetPhysState (void) const
 This function indicates whether the array of physical values $\boldsymbol{u}_l$ (implemented as m_phys) is filled or not.
NekDouble PhysIntegral (void)
 This function integrates a function $f(\boldsymbol{x})$ over the domain consisting of all the elements of the expansion.
NekDouble PhysIntegral (const Array< OneD, const NekDouble > &inarray)
 This function integrates a function $f(\boldsymbol{x})$ over the domain consisting of all the elements of the expansion.
void IProductWRTBase_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function calculates the inner product of a function $f(\boldsymbol{x})$ with respect to all {local} expansion modes $\phi_n^e(\boldsymbol{x})$.
void IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
void IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function calculates the inner product of a function $f(\boldsymbol{x})$ with respect to the derivative (in direction.
void FwdTrans_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function elementally evaluates the forward transformation of a function $u(\boldsymbol{x})$ onto the global spectral/hp expansion.
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
void MultiplyByElmtInvMass (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function elementally mulplies the coefficient space of Sin my the elemental inverse of the mass matrix.
void MultiplyByInvMassMatrix (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
void SmoothField (Array< OneD, NekDouble > &field)
 Smooth a field across elements.
void HelmSolve (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const FlagList &flags, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff=StdRegions::NullVarCoeffMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 Solve helmholtz problem.
void LinearAdvectionDiffusionReactionSolve (const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 Solve Advection Diffusion Reaction.
void LinearAdvectionReactionSolve (const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 Solve Advection Diffusion Reaction.
void FwdTrans_BndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void BwdTrans_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function elementally evaluates the backward transformation of the global spectral/hp element expansion.
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
void GetCoords (Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1=NullNekDouble1DArray, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
 This function calculates the coordinates of all the elemental quadrature points $\boldsymbol{x}_i$.
void 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 ApplyGeomInfo ()
 Apply geometry information to each expansion.
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)
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.
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 (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.
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
boost::shared_ptr< ExpList > & UpdateBndCondExpansion (int i)
void Upwind (const Array< OneD, const Array< OneD, NekDouble > > &Vec, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
void Upwind (const Array< OneD, const NekDouble > &Vn, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
boost::shared_ptr< ExpList > & GetTrace ()
boost::shared_ptr
< AssemblyMapDG > & 
GetTraceMap (void)
const Array< OneD, const int > & GetTraceBndMap (void)
void GetNormals (Array< OneD, Array< OneD, NekDouble > > &normals)
void AddTraceIntegral (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
void AddTraceIntegral (const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
void AddFwdBwdTraceIntegral (const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
void GetFwdBwdTracePhys (Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
void GetFwdBwdTracePhys (const Array< OneD, const NekDouble > &field, Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
void ExtractTracePhys (Array< OneD, NekDouble > &outarray)
void ExtractTracePhys (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
const Array< OneD, const
SpatialDomains::BoundaryConditionShPtr > & 
GetBndConditions ()
Array< OneD,
SpatialDomains::BoundaryConditionShPtr > & 
UpdateBndConditions ()
void EvaluateBoundaryConditions (const NekDouble time=0.0, const std::string varName="", const NekDouble=NekConstants::kNekUnsetDouble, const NekDouble=NekConstants::kNekUnsetDouble)
void GeneralMatrixOp (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
 This function calculates the result of the multiplication of a matrix of type specified by mkey with a vector given by inarray.
void GeneralMatrixOp_IterPerExp (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void SetUpPhysNormals ()
void GetBoundaryToElmtMap (Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
void GeneralGetFieldDefinitions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef, int NumHomoDir=0, Array< OneD, LibUtilities::BasisSharedPtr > &HomoBasis=LibUtilities::NullBasisSharedPtr1DArray, std::vector< NekDouble > &HomoLen=LibUtilities::NullNekDoubleVector, std::vector< unsigned int > &HomoZIDs=LibUtilities::NullUnsignedIntVector, std::vector< unsigned int > &HomoYIDs=LibUtilities::NullUnsignedIntVector)
const
NekOptimize::GlobalOptParamSharedPtr
GetGlobalOptParam (void)
map< int, RobinBCInfoSharedPtrGetRobinBCInfo ()
void GetPeriodicEntities (PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces=NullPeriodicMap)
std::vector
< LibUtilities::FieldDefinitionsSharedPtr
GetFieldDefinitions ()
void GetFieldDefinitions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef)
void AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
 Append the element data listed in elements fielddef->m_ElementIDs onto fielddata.
void AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, Array< OneD, NekDouble > &coeffs)
 Append the data in coeffs listed in elements fielddef->m_ElementIDs onto fielddata.
void ExtractElmtDataToCoeffs (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
 Extract the data in fielddata into the coeffs using the basic ExpList Elemental expansions rather than planes in homogeneous case.
void ExtractCoeffsToCoeffs (const boost::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
 Extract the data from fromField using fromExpList the coeffs using the basic ExpList Elemental expansions rather than planes in homogeneous case.
void ExtractDataToCoeffs (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
 Extract the data in fielddata into the coeffs.
boost::shared_ptr< ExpListGetSharedThisPtr ()
 Returns a shared pointer to the current object.
boost::shared_ptr
< LibUtilities::SessionReader
GetSession ()
 Returns the session object.
boost::shared_ptr
< LibUtilities::Comm
GetComm ()
 Returns the comm object.
SpatialDomains::MeshGraphSharedPtr GetGraph ()
LibUtilities::BasisSharedPtr GetHomogeneousBasis (void)
boost::shared_ptr< ExpList > & GetPlane (int n)

Protected Member Functions

virtual void v_BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
 Performs the backward transformation of the spectral/hp element expansion.
virtual void v_IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
 Calculates the inner product of a function $f(\boldsymbol{x})$ with respect to all global expansion modes $\phi_n^e(\boldsymbol{x})$.
virtual void v_FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
- Protected Member Functions inherited from Nektar::MultiRegions::DisContField3D
void SetUpDG (const std::string="DefaultVar")
 Set up all DG member variables and maps.
bool SameTypeOfBoundaryConditions (const DisContField3D &In)
void GenerateBoundaryConditionExpansion (const SpatialDomains::MeshGraphSharedPtr &graph3D, const SpatialDomains::BoundaryConditions &bcs, const std::string &variable)
void FindPeriodicFaces (const SpatialDomains::BoundaryConditions &bcs, const std::string &variable)
 Determine the periodic faces, edges and vertices for the given graph.
bool IsLeftAdjacentFace (const int n, const int e)
virtual void v_GetFwdBwdTracePhys (Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
 This method extracts the "forward" and "backward" trace data from the array field and puts the data into output vectors Fwd and 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_AddTraceIntegral (const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 Add trace contributions into elemental coefficient spaces.
virtual void v_AddFwdBwdTraceIntegral (const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
 Add trace contributions into elemental coefficient spaces.
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_GeneralMatrixOp (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
 Calculates the result of the multiplication of a global matrix of type specified by mkey with a vector given by inarray.
virtual void v_GetBoundaryToElmtMap (Array< OneD, int > &ElmtID, Array< OneD, int > &FaceID)
 Set up a list of elemeent IDs and edge IDs that link to the boundary conditions.
virtual void v_GetPeriodicEntities (PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces)
virtual ExpListSharedPtrv_GetTrace ()
virtual AssemblyMapDGSharedPtrv_GetTraceMap ()
virtual const Array< OneD,
const
MultiRegions::ExpListSharedPtr > & 
v_GetBndCondExpansions ()
virtual const Array< OneD,
const
SpatialDomains::BoundaryConditionShPtr > & 
v_GetBndConditions ()
virtual
MultiRegions::ExpListSharedPtr
v_UpdateBndCondExpansion (int i)
virtual Array< OneD,
SpatialDomains::BoundaryConditionShPtr > & 
v_UpdateBndConditions ()
virtual void v_EvaluateBoundaryConditions (const NekDouble time=0.0, const std::string varName="", const NekDouble x2_in=NekConstants::kNekUnsetDouble, const NekDouble x3_in=NekConstants::kNekUnsetDouble)
 This function evaluates the boundary conditions at a certain time-level.
virtual map< int,
RobinBCInfoSharedPtr
v_GetRobinBCInfo ()
- Protected Member Functions inherited from Nektar::MultiRegions::ExpList3D
virtual void v_SetUpPhysNormals ()
 Set up the normals on each expansion.
- 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 int v_GetNumElmts (void)
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 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_MultiplyByInvMassMatrix (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
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_BwdTrans_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_FwdTrans_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_SmoothField (Array< OneD, NekDouble > &field)
virtual void v_IProductWRTBase_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_GetCoords (Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
virtual void v_PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
virtual void v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
virtual void v_PhysDeriv (Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
virtual void v_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_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_ReadGlobalOptimizationParameters ()
virtual std::vector
< LibUtilities::FieldDefinitionsSharedPtr
v_GetFieldDefinitions (void)
virtual void v_GetFieldDefinitions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef)
virtual void v_AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
virtual void v_AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, Array< OneD, NekDouble > &coeffs)
virtual void v_ExtractDataToCoeffs (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
 Extract data from raw field data into expansion list.
virtual void v_ExtractCoeffsToCoeffs (const boost::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
virtual void v_WriteTecplotHeader (std::ostream &outfile, std::string var="")
virtual void v_WriteTecplotZone (std::ostream &outfile, int expansion)
virtual void v_WriteTecplotField (std::ostream &outfile, int expansion)
virtual void v_WriteTecplotConnectivity (std::ostream &outfile, int expansion)
virtual void v_WriteVtkPieceHeader (std::ostream &outfile, int expansion)
virtual void v_WriteVtkPieceData (std::ostream &outfile, int expansion, std::string var)
virtual NekDouble v_L2 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
virtual NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray)
virtual Array< OneD, const
NekDouble
v_HomogeneousEnergy (void)
virtual
LibUtilities::TranspositionSharedPtr 
v_GetTransposition (void)
virtual NekDouble v_GetHomoLen (void)
virtual Array< OneD, const
unsigned int > 
v_GetZIDs (void)
virtual Array< OneD, const
unsigned int > 
v_GetYIDs (void)
virtual void v_PhysInterp1DScaled (const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_PhysGalerkinProjection1DScaled (const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void ExtractFileBCs (const std::string &fileName, const std::string &varName, const boost::shared_ptr< ExpList > locExpList)

Protected Attributes

AssemblyMapCGSharedPtr m_locToGloMap
GlobalMatrixMapShPtr m_globalMat
 (A shared pointer to) a list which collects all the global matrices being assembled, such that they should be constructed only once.
LibUtilities::NekManager
< GlobalLinSysKey,
GlobalLinSys
m_globalLinSysManager
 (A shared pointer to) a list which collects all the global linear system being assembled, such that they should be constructed only once.
- Protected Attributes inherited from Nektar::MultiRegions::DisContField3D
Array< OneD,
MultiRegions::ExpListSharedPtr
m_bndCondExpansions
 An object which contains the discretised boundary conditions.
Array< OneD,
SpatialDomains::BoundaryConditionShPtr
m_bndConditions
 An array which contains the information about the boundary condition on the different boundary regions.
GlobalLinSysMapShPtr m_globalBndMat
ExpListSharedPtr m_trace
AssemblyMapDGSharedPtr m_traceMap
std::set< int > m_boundaryFaces
 A set storing the global IDs of any boundary faces.
PeriodicMap m_periodicFaces
 A map which identifies pairs of periodic faces.
PeriodicMap m_periodicEdges
 A map which identifies groups of periodic edges.
PeriodicMap m_periodicVerts
 A map which identifies groups of periodic vertices.
vector< bool > m_leftAdjacentFaces
vector< int > m_periodicFwdCopy
 A vector indicating degress of freedom which need to be copied from forwards to backwards space in case of a periodic boundary condition.
vector< int > m_periodicBwdCopy

Private Member Functions

GlobalLinSysSharedPtr GetGlobalLinSys (const GlobalLinSysKey &mkey)
GlobalLinSysSharedPtr GenGlobalLinSys (const GlobalLinSysKey &mkey)
GlobalMatrixSharedPtr GetGlobalMatrix (const GlobalMatrixKey &mkey)
 Returns the global matrix specified by mkey.
void GlobalSolve (const GlobalLinSysKey &key, const Array< OneD, const NekDouble > &rhs, Array< OneD, NekDouble > &inout, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
virtual void v_ImposeDirichletConditions (Array< OneD, NekDouble > &outarray)
 Impose the Dirichlet Boundary Conditions on outarray.
virtual void v_FillBndCondFromField ()
virtual void v_LocalToGlobal (void)
virtual void v_GlobalToLocal (void)
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_GeneralMatrixOp (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)

Additional Inherited Members

- 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)

Detailed Description

Definition at line 51 of file ContField3D.h.

Constructor & Destructor Documentation

Nektar::MultiRegions::ContField3D::ContField3D ( )

Definition at line 45 of file ContField3D.cpp.

:
boost::bind(&ContField3D::GenGlobalLinSys, this, _1),
std::string("GlobalLinSys"))
{
}
Nektar::MultiRegions::ContField3D::ContField3D ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr graph3D,
const std::string &  variable = "DefaultVar",
const bool  CheckIfSingularSystem = false 
)

Construct a global continuous field.

Given a mesh graph2D, containing information about the domain and the spectral/hp element expansion, this constructor fills the list of local expansions m_exp with the proper expansions, calculates the total number of quadrature points $\boldsymbol{x}_i$ and local expansion coefficients $\hat{u}^e_n$ and allocates memory for the arrays m_coeffs and m_phys. Furthermore, it constructs the mapping array (contained in m_locToGloMap) for the transformation between local elemental level and global level, it calculates the total number global expansion coefficients $\hat{u}_n$ and allocates memory for the array m_coeffs. The constructor also discretises the boundary conditions, specified by the argument bcs, by expressing them in terms of the coefficient of the expansion on the boundary.

Parameters
pSessionSession information.
graph3DA mesh, containing information about the domain and the spectral/hp element expansion.
variableThe variable associated with this field.

Definition at line 76 of file ContField3D.cpp.

References Nektar::MultiRegions::DisContField3D::m_bndCondExpansions, Nektar::MultiRegions::DisContField3D::m_bndConditions, m_locToGloMap, Nektar::MultiRegions::ExpList::m_ncoeffs, Nektar::MultiRegions::DisContField3D::m_periodicEdges, Nektar::MultiRegions::DisContField3D::m_periodicFaces, Nektar::MultiRegions::DisContField3D::m_periodicVerts, and Nektar::MultiRegions::ExpList::m_session.

:
DisContField3D(pSession,graph3D,variable,false),
boost::bind(&ContField3D::GenGlobalLinSys, this, _1),
std::string("GlobalLinSys"))
{
CheckIfSingularSystem, variable,
if (m_session->DefinesCmdLineArgument("verbose"))
{
m_locToGloMap->PrintStats(std::cout, variable);
}
}
Nektar::MultiRegions::ContField3D::ContField3D ( const ContField3D In,
const SpatialDomains::MeshGraphSharedPtr graph3D,
const std::string &  variable,
const bool  CheckIfSingularSystem = false 
)

Construct a global continuous field with solution type based on another field but using a separate input mesh and boundary conditions.

Given a mesh graph3D, containing information about the domain and the spectral/hp element expansion, this constructor fills the list of local expansions m_exp with the proper expansions, calculates the total number of quadrature points $\boldsymbol{x}_i$ and local expansion coefficients $\hat{u}^e_n$ and allocates memory for the arrays m_coeffs and m_phys. Furthermore, it constructs the mapping array (contained in m_locToGloMap) for the transformation between local elemental level and global level, it calculates the total number global expansion coefficients $\hat{u}_n$ and allocates memory for the array m_coeffs. The constructor also discretises the boundary conditions, specified by the argument bcs, by expressing them in terms of the coefficient of the expansion on the boundary.

Parameters
InExisting ContField2D object used to provide the local to global mapping information and global solution type.
graph3DA mesh, containing information about the domain and the spectral/hp element expansion.
bcsThe boundary conditions.
bc_loc

Definition at line 121 of file ContField3D.cpp.

References Nektar::MultiRegions::DisContField3D::m_bndCondExpansions, Nektar::MultiRegions::DisContField3D::m_bndConditions, m_locToGloMap, Nektar::MultiRegions::ExpList::m_ncoeffs, Nektar::MultiRegions::DisContField3D::m_periodicEdges, Nektar::MultiRegions::DisContField3D::m_periodicFaces, Nektar::MultiRegions::DisContField3D::m_periodicVerts, Nektar::MultiRegions::ExpList::m_session, and Nektar::MultiRegions::DisContField3D::SameTypeOfBoundaryConditions().

:
DisContField3D(In,graph3D,variable,false),
std::string("GlobalLinSys"))
{
if(!SameTypeOfBoundaryConditions(In) || CheckIfSingularSystem)
{
CheckIfSingularSystem, variable,
if (m_session->DefinesCmdLineArgument("verbose"))
{
m_locToGloMap->PrintStats(std::cout, variable);
}
}
else
{
m_locToGloMap = In.m_locToGloMap;
}
}
Nektar::MultiRegions::ContField3D::ContField3D ( const ContField3D In)

Definition at line 151 of file ContField3D.cpp.

:
m_locToGloMap(In.m_locToGloMap),
m_globalMat(In.m_globalMat),
m_globalLinSysManager(In.m_globalLinSysManager)
{
}
Nektar::MultiRegions::ContField3D::~ContField3D ( )
virtual

Definition at line 160 of file ContField3D.cpp.

{
}

Member Function Documentation

void Nektar::MultiRegions::ContField3D::Assemble ( )
inline

Definition at line 215 of file ContField3D.h.

References Nektar::MultiRegions::ExpList::m_coeffs, and m_locToGloMap.

Referenced by v_GeneralMatrixOp(), v_IProductWRTBase(), and v_MultiplyByInvMassMatrix().

{
}
void Nektar::MultiRegions::ContField3D::Assemble ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
inline

Definition at line 220 of file ContField3D.h.

References m_locToGloMap.

{
m_locToGloMap->Assemble(inarray, outarray);
}
void Nektar::MultiRegions::ContField3D::GenerateDirBndCondForcing ( const GlobalLinSysKey key,
Array< OneD, NekDouble > &  inout,
Array< OneD, NekDouble > &  outarray 
)

Definition at line 323 of file ContField3D.cpp.

References Nektar::SpatialDomains::eDirichlet, Nektar::MultiRegions::eGlobal, Nektar::MultiRegions::ExpList::GeneralMatrixOp(), Nektar::MultiRegions::ExpList::GetNcoeffs(), Nektar::MultiRegions::DisContField3D::m_bndCondExpansions, Nektar::MultiRegions::DisContField3D::m_bndConditions, m_locToGloMap, and sign.

{
int bndcnt=0;
const Array<OneD,const int>& map = m_locToGloMap->GetBndCondCoeffsToGlobalCoeffsMap();
for(int i = 0; i < m_bndCondExpansions.num_elements(); ++i)
{
if(m_bndConditions[i]->GetBoundaryConditionType() == SpatialDomains::eDirichlet)
{
const Array<OneD,const NekDouble>& coeffs = m_bndCondExpansions[i]->GetCoeffs();
for(int j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); ++j)
{
sign = m_locToGloMap->GetBndCondCoeffsToGlobalCoeffsSign(bndcnt);
inout[map[bndcnt++]] = sign * coeffs[j];
}
}
else
{
bndcnt += m_bndCondExpansions[i]->GetNcoeffs();
}
}
GeneralMatrixOp(key,inout,outarray,eGlobal);
}
GlobalLinSysSharedPtr Nektar::MultiRegions::ContField3D::GenGlobalLinSys ( const GlobalLinSysKey mkey)
private

Definition at line 381 of file ContField3D.cpp.

References ASSERTL1, Nektar::MultiRegions::GlobalMatrixKey::LocToGloMapIsDefined(), and m_locToGloMap.

{
ASSERTL1(mkey.LocToGloMapIsDefined(),
"To use method must have a AssemblyMap "
"attached to key");
}
const Array< OneD, const MultiRegions::ExpListSharedPtr > & Nektar::MultiRegions::ContField3D::GetBndCondExp ( )
inline

This function return the boundary conditions expansion.

Definition at line 193 of file ContField3D.h.

References Nektar::MultiRegions::DisContField3D::m_bndCondExpansions.

{
}
const Array<OneD,const MultiRegions::ExpListSharedPtr>& Nektar::MultiRegions::ContField3D::GetBndCondExpansions ( )
inline

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 75 of file ContField3D.h.

References Nektar::MultiRegions::DisContField3D::m_bndCondExpansions.

{
}
GlobalLinSysSharedPtr Nektar::MultiRegions::ContField3D::GetGlobalLinSys ( const GlobalLinSysKey mkey)
private

Definition at line 375 of file ContField3D.cpp.

References m_globalLinSysManager.

Referenced by GlobalSolve().

{
return m_globalLinSysManager[mkey];
}
GlobalMatrixSharedPtr Nektar::MultiRegions::ContField3D::GetGlobalMatrix ( const GlobalMatrixKey mkey)
private

Returns the global matrix specified by mkey.

Returns the global matrix associated with the given GlobalMatrixKey. If the global matrix has not yet been constructed on this field, it is first constructed using GenGlobalMatrix().

Parameters
mkeyGlobal matrix key.
Returns
Assocated global matrix.

Definition at line 397 of file ContField3D.cpp.

References ASSERTL1, Nektar::MultiRegions::ExpList::GenGlobalMatrix(), Nektar::iterator, Nektar::MultiRegions::GlobalMatrixKey::LocToGloMapIsDefined(), m_globalMat, and m_locToGloMap.

Referenced by v_BwdTrans(), v_GeneralMatrixOp(), and v_IProductWRTBase().

{
ASSERTL1(mkey.LocToGloMapIsDefined(),
"To use method must have a AssemblyMap "
"attached to key");
GlobalMatrixMap::iterator matrixIter = m_globalMat->find(mkey);
if(matrixIter == m_globalMat->end())
{
glo_matrix = GenGlobalMatrix(mkey,m_locToGloMap);
(*m_globalMat)[mkey] = glo_matrix;
}
else
{
glo_matrix = matrixIter->second;
}
return glo_matrix;
}
int Nektar::MultiRegions::ContField3D::GetGlobalMatrixNnz ( const GlobalMatrixKey gkey)

Definition at line 602 of file ContField3D.cpp.

References ASSERTL1, Nektar::iterator, Nektar::MultiRegions::GlobalMatrixKey::LocToGloMapIsDefined(), and m_globalMat.

{
ASSERTL1(gkey.LocToGloMapIsDefined(),
"To use method must have a AssemblyMap "
"attached to key");
GlobalMatrixMap::iterator matrixIter = m_globalMat->find(gkey);
if(matrixIter == m_globalMat->end())
{
return 0;
}
else
{
return matrixIter->second->GetNumNonZeroEntries();
}
return 0;
}
const AssemblyMapCGSharedPtr & Nektar::MultiRegions::ContField3D::GetLocalToGlobalMap ( ) const
inline

Definition at line 228 of file ContField3D.h.

References m_locToGloMap.

{
return m_locToGloMap;
}
void Nektar::MultiRegions::ContField3D::GlobalSolve ( const GlobalLinSysKey key,
const Array< OneD, const NekDouble > &  rhs,
Array< OneD, NekDouble > &  inout,
const Array< OneD, const NekDouble > &  dirForcing = NullNekDouble1DArray 
)
private

Definition at line 353 of file ContField3D.cpp.

References GetGlobalLinSys(), m_locToGloMap, and v_ImposeDirichletConditions().

Referenced by v_FwdTrans(), v_HelmSolve(), and v_MultiplyByInvMassMatrix().

{
int NumDirBcs = m_locToGloMap->GetNumGlobalDirBndCoeffs();
int contNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
// STEP 1: SET THE DIRICHLET DOFS TO THE RIGHT VALUE
// IN THE SOLUTION ARRAY
// STEP 2: CALCULATE THE HOMOGENEOUS COEFFICIENTS
if(contNcoeffs - NumDirBcs > 0)
{
LinSys->Solve(rhs,inout,m_locToGloMap,dirForcing);
}
}
void Nektar::MultiRegions::ContField3D::GlobalToLocal ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
inline

Definition at line 200 of file ContField3D.h.

References m_locToGloMap.

{
m_locToGloMap->GlobalToLocal(inarray, outarray);
}
void Nektar::MultiRegions::ContField3D::LocalToGlobal ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
inline

Definition at line 208 of file ContField3D.h.

References m_locToGloMap.

{
m_locToGloMap->LocalToGlobal(inarray, outarray);
}
void Nektar::MultiRegions::ContField3D::v_BwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
CoeffState  coeffstate = eLocal 
)
protectedvirtual

Performs the backward transformation of the spectral/hp element expansion.

Given the coefficients of an expansion, this function evaluates the spectral/hp expansion $u^{\delta}(\boldsymbol{x})$ at the quadrature points $\boldsymbol{x}_i$. This operation is evaluated locally by the function ExpList::BwdTrans.

The coefficients of the expansion should be contained in the variable #inarray of the ExpList object In. The resulting physical values at the quadrature points $u^{\delta}(\boldsymbol{x}_i)$ are stored in the array #outarray.

Parameters
InAn ExpList, containing the local coefficients $\hat{u}_n^e$ in its array m_coeffs.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 179 of file ContField3D.cpp.

References Nektar::MultiRegions::ExpList::BwdTrans_IterPerExp(), Nektar::StdRegions::eBwdTrans, Nektar::MultiRegions::eGlobal, GetGlobalMatrix(), Nektar::MultiRegions::ExpList::GlobalToLocal(), Nektar::MultiRegions::ExpList::m_globalOptParam, m_locToGloMap, and Nektar::MultiRegions::ExpList::m_ncoeffs.

{
if(coeffstate == eGlobal)
{
bool doGlobalOp = m_globalOptParam->DoGlobalMatOp(
if(doGlobalOp)
{
GlobalMatrixKey gkey(StdRegions::eBwdTrans,m_locToGloMap);
mat->Multiply(inarray,outarray);
}
else
{
Array<OneD, NekDouble> wsp(m_ncoeffs);
GlobalToLocal(inarray,wsp);
BwdTrans_IterPerExp(wsp,outarray);
}
}
else
{
BwdTrans_IterPerExp(inarray,outarray);
}
}
void Nektar::MultiRegions::ContField3D::v_FillBndCondFromField ( void  )
privatevirtual

Definition at line 477 of file ContField3D.cpp.

References Nektar::MultiRegions::ExpList::GetNcoeffs(), Nektar::MultiRegions::ExpList::LocalToGlobal(), Nektar::MultiRegions::DisContField3D::m_bndCondExpansions, Nektar::MultiRegions::ExpList::m_coeffs, m_locToGloMap, and sign.

{
int bndcnt = 0;
const Array<OneD,const int> &bndMap =
m_locToGloMap->GetBndCondCoeffsToGlobalCoeffsMap();
Array<OneD, NekDouble> tmp(m_locToGloMap->GetNumGlobalCoeffs());
// Now fill in all other Dirichlet coefficients.
for(int i = 0; i < m_bndCondExpansions.num_elements(); ++i)
{
Array<OneD, NekDouble>& coeffs = m_bndCondExpansions[i]->UpdateCoeffs();
for(int j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); ++j)
{
sign = m_locToGloMap->GetBndCondCoeffsToGlobalCoeffsSign(bndcnt);
coeffs[j] = sign * tmp[bndMap[bndcnt++]];
}
}
}
void Nektar::MultiRegions::ContField3D::v_FwdTrans ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
CoeffState  coeffstate 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 255 of file ContField3D.cpp.

References Nektar::MultiRegions::eGlobal, Nektar::StdRegions::eMass, GlobalSolve(), Nektar::MultiRegions::ExpList::GlobalToLocal(), Nektar::MultiRegions::ExpList::IProductWRTBase(), and m_locToGloMap.

{
// Inner product of forcing
int contNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
Array<OneD,NekDouble> wsp(contNcoeffs);
IProductWRTBase(inarray,wsp,eGlobal);
// Solve the system
GlobalLinSysKey key(StdRegions::eMass, m_locToGloMap);
if(coeffstate == eGlobal)
{
GlobalSolve(key,wsp,outarray);
}
else
{
Array<OneD,NekDouble> tmp(contNcoeffs,0.0);
GlobalSolve(key,wsp,tmp);
GlobalToLocal(tmp,outarray);
}
}
void Nektar::MultiRegions::ContField3D::v_GeneralMatrixOp ( const GlobalMatrixKey gkey,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
CoeffState  coeffstate 
)
privatevirtual

Definition at line 571 of file ContField3D.cpp.

References Assemble(), Nektar::MultiRegions::eGlobal, Nektar::MultiRegions::ExpList::GeneralMatrixOp_IterPerExp(), GetGlobalMatrix(), Nektar::MultiRegions::GlobalMatrixKey::GetMatrixType(), Nektar::MultiRegions::ExpList::GlobalToLocal(), Nektar::MultiRegions::ExpList::m_globalOptParam, m_locToGloMap, and Nektar::MultiRegions::ExpList::m_ncoeffs.

{
if(coeffstate == eGlobal)
{
bool doGlobalOp = m_globalOptParam->DoGlobalMatOp(gkey.GetMatrixType());
if(doGlobalOp)
{
mat->Multiply(inarray,outarray);
m_locToGloMap->UniversalAssemble(outarray);
}
else
{
Array<OneD,NekDouble> tmp1(2*m_ncoeffs);
Array<OneD,NekDouble> tmp2(tmp1+m_ncoeffs);
GlobalToLocal(inarray,tmp1);
GeneralMatrixOp_IterPerExp(gkey,tmp1,tmp2);
Assemble(tmp2,outarray);
}
}
else
{
GeneralMatrixOp_IterPerExp(gkey,inarray,outarray);
}
}
void Nektar::MultiRegions::ContField3D::v_GlobalToLocal ( void  )
privatevirtual

Definition at line 505 of file ContField3D.cpp.

References Nektar::MultiRegions::ExpList::m_coeffs, and m_locToGloMap.

{
m_locToGloMap->GlobalToLocal(m_coeffs, m_coeffs);
}
void Nektar::MultiRegions::ContField3D::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 512 of file ContField3D.cpp.

References Nektar::SpatialDomains::eDirichlet, Nektar::MultiRegions::eGlobal, Nektar::StdRegions::eHelmholtz, Nektar::eUseGlobal, Nektar::MultiRegions::ExpList::GetNcoeffs(), GlobalSolve(), Nektar::MultiRegions::ExpList::GlobalToLocal(), Nektar::MultiRegions::ExpList::IProductWRTBase(), Nektar::FlagList::isSet(), Nektar::MultiRegions::ExpList::LocalToGlobal(), Nektar::MultiRegions::DisContField3D::m_bndCondExpansions, Nektar::MultiRegions::DisContField3D::m_bndConditions, m_locToGloMap, Vmath::Neg(), sign, and Vmath::Vadd().

{
// Inner product of forcing
int contNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
Array<OneD,NekDouble> wsp(contNcoeffs);
IProductWRTBase(inarray,wsp,eGlobal);
// Note -1.0 term necessary to invert forcing function to
// be consistent with matrix definition
Vmath::Neg(contNcoeffs, wsp, 1);
// Forcing function with weak boundary conditions
int i,j;
int bndcnt = 0;
Array<OneD, NekDouble> gamma(contNcoeffs, 0.0);
for(i = 0; i < m_bndCondExpansions.num_elements(); ++i)
{
if(m_bndConditions[i]->GetBoundaryConditionType() != SpatialDomains::eDirichlet)
{
for(j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); j++)
{
sign = m_locToGloMap->GetBndCondCoeffsToGlobalCoeffsSign(bndcnt);
gamma[m_locToGloMap->GetBndCondCoeffsToGlobalCoeffsMap(bndcnt++)] +=
sign * (m_bndCondExpansions[i]->GetCoeffs())[j];
}
}
else
{
bndcnt += m_bndCondExpansions[i]->GetNcoeffs();
}
}
m_locToGloMap->UniversalAssemble(gamma);
// Add weak boundary conditions to forcing
Vmath::Vadd(contNcoeffs, wsp, 1, gamma, 1, wsp, 1);
// Solve the system
GlobalLinSysKey key(StdRegions::eHelmholtz, m_locToGloMap, factors,varcoeff);
if(flags.isSet(eUseGlobal))
{
GlobalSolve(key,wsp,outarray,dirForcing);
}
else
{
Array<OneD,NekDouble> tmp(contNcoeffs);
LocalToGlobal(outarray,tmp);
GlobalSolve(key,wsp,tmp,dirForcing);
GlobalToLocal(tmp,outarray);
}
}
void Nektar::MultiRegions::ContField3D::v_ImposeDirichletConditions ( Array< OneD, NekDouble > &  outarray)
privatevirtual

Impose the Dirichlet Boundary Conditions on outarray.

Definition at line 420 of file ContField3D.cpp.

References Nektar::SpatialDomains::eDirichlet, Nektar::MultiRegions::ExpList::GetNcoeffs(), Nektar::iterator, Nektar::MultiRegions::DisContField3D::m_bndCondExpansions, Nektar::MultiRegions::DisContField3D::m_bndConditions, m_locToGloMap, sign, and Vmath::Vcopy().

Referenced by GlobalSolve().

{
int i,j;
int bndcnt = 0;
int nDir = m_locToGloMap->GetNumGlobalDirBndCoeffs();
const Array<OneD,const int> &bndMap =
m_locToGloMap->GetBndCondCoeffsToGlobalCoeffsMap();
Array<OneD, NekDouble> tmp(
m_locToGloMap->GetNumGlobalBndCoeffs(), 0.0);
// Fill in Dirichlet coefficients that are to be sent to other
// processors. This code block uses a tuple<int,int.NekDouble> which
// stores the local id of coefficent the global id of the data
// location and the inverse of the values of the data (arising from
// periodic boundary conditiosn)
map<int, vector<ExtraDirDof> > &extraDirDofs =
m_locToGloMap->GetExtraDirDofs();
map<int, vector<ExtraDirDof> >::iterator it;
for (it = extraDirDofs.begin(); it != extraDirDofs.end(); ++it)
{
for (i = 0; i < it->second.size(); ++i)
{
tmp[it->second.at(i).get<1>()] =
m_bndCondExpansions[it->first]->GetCoeffs()[
it->second.at(i).get<0>()]*it->second.at(i).get<2>();
}
}
m_locToGloMap->UniversalAssembleBnd(tmp);
// Now fill in all other Dirichlet coefficients.
for(i = 0; i < m_bndCondExpansions.num_elements(); ++i)
{
if(m_bndConditions[i]->GetBoundaryConditionType() ==
{
const Array<OneD,const NekDouble>& coeffs =
m_bndCondExpansions[i]->GetCoeffs();
for(j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); ++j)
{
sign = m_locToGloMap->GetBndCondCoeffsToGlobalCoeffsSign(
bndcnt);
tmp[bndMap[bndcnt++]] = sign * coeffs[j];
}
}
else
{
bndcnt += m_bndCondExpansions[i]->GetNcoeffs();
}
}
Vmath::Vcopy(nDir, tmp, 1, outarray, 1);
}
void Nektar::MultiRegions::ContField3D::v_IProductWRTBase ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
CoeffState  coeffstate = eLocal 
)
protectedvirtual

Calculates the inner product of a function $f(\boldsymbol{x})$ with respect to all global expansion modes $\phi_n^e(\boldsymbol{x})$.

The operation is evaluated locally (i.e. with respect to all local expansion modes) by the function ExpList::IProductWRTBase. The inner product with respect to the global expansion modes is than obtained by a global assembly operation.

The values of the function $f(\boldsymbol{x})$ evaluated at the quadrature points $\boldsymbol{x}_i$ should be contained in the variable m_phys of the ExpList object in. The result is stored in the array m_coeffs.

Parameters
InAn ExpList, containing the discrete evaluation of $f(\boldsymbol{x})$ at the quadrature points in its array m_phys.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 224 of file ContField3D.cpp.

References Assemble(), Nektar::MultiRegions::eGlobal, Nektar::StdRegions::eIProductWRTBase, GetGlobalMatrix(), Nektar::MultiRegions::ExpList::IProductWRTBase_IterPerExp(), Nektar::MultiRegions::ExpList::m_globalOptParam, m_locToGloMap, and Nektar::MultiRegions::ExpList::m_ncoeffs.

{
if(coeffstate == eGlobal)
{
bool doGlobalOp = m_globalOptParam->DoGlobalMatOp(
if(doGlobalOp)
{
GlobalMatrixKey gkey(StdRegions::eIProductWRTBase,
mat->Multiply(inarray,outarray);
m_locToGloMap->UniversalAssemble(outarray);
}
else
{
Array<OneD, NekDouble> wsp(m_ncoeffs);
Assemble(wsp,outarray);
}
}
else
{
IProductWRTBase_IterPerExp(inarray,outarray);
}
}
void Nektar::MultiRegions::ContField3D::v_LocalToGlobal ( void  )
privatevirtual

Definition at line 500 of file ContField3D.cpp.

References Nektar::MultiRegions::ExpList::m_coeffs, and m_locToGloMap.

{
m_locToGloMap->LocalToGlobal(m_coeffs, m_coeffs);
}
void Nektar::MultiRegions::ContField3D::v_MultiplyByInvMassMatrix ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
CoeffState  coeffstate 
)
privatevirtual

Definition at line 280 of file ContField3D.cpp.

References Assemble(), Nektar::MultiRegions::eGlobal, Nektar::StdRegions::eMass, GlobalSolve(), Nektar::MultiRegions::ExpList::GlobalToLocal(), m_locToGloMap, and Vmath::Vcopy().

{
int contNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
GlobalLinSysKey key(StdRegions::eMass, m_locToGloMap);
if(coeffstate == eGlobal)
{
if(inarray.data() == outarray.data())
{
Array<OneD, NekDouble> tmp(contNcoeffs,0.0);
Vmath::Vcopy(contNcoeffs,inarray,1,tmp,1);
GlobalSolve(key,tmp,outarray);
}
else
{
GlobalSolve(key,inarray,outarray);
}
}
else
{
Array<OneD, NekDouble> globaltmp(contNcoeffs,0.0);
if(inarray.data() == outarray.data())
{
Array<OneD,NekDouble> tmp(inarray.num_elements());
Vmath::Vcopy(inarray.num_elements(),inarray,1,tmp,1);
Assemble(tmp,outarray);
}
else
{
Assemble(inarray,outarray);
}
GlobalSolve(key,outarray,globaltmp);
GlobalToLocal(globaltmp,outarray);
}
}

Member Data Documentation

LibUtilities::NekManager<GlobalLinSysKey, GlobalLinSys> Nektar::MultiRegions::ContField3D::m_globalLinSysManager
protected

(A shared pointer to) a list which collects all the global linear system being assembled, such that they should be constructed only once.

Definition at line 120 of file ContField3D.h.

Referenced by GetGlobalLinSys().

GlobalMatrixMapShPtr Nektar::MultiRegions::ContField3D::m_globalMat
protected

(A shared pointer to) a list which collects all the global matrices being assembled, such that they should be constructed only once.

Definition at line 115 of file ContField3D.h.

Referenced by GetGlobalMatrix(), and GetGlobalMatrixNnz().

AssemblyMapCGSharedPtr Nektar::MultiRegions::ContField3D::m_locToGloMap
protected

Definition at line 110 of file ContField3D.h.

Referenced by Assemble(), ContField3D(), GenerateDirBndCondForcing(), GenGlobalLinSys(), GetGlobalMatrix(), GetLocalToGlobalMap(), GlobalSolve(), GlobalToLocal(), LocalToGlobal(), v_BwdTrans(), v_FillBndCondFromField(), v_FwdTrans(), v_GeneralMatrixOp(), v_GlobalToLocal(), v_HelmSolve(), v_ImposeDirichletConditions(), v_IProductWRTBase(), v_LocalToGlobal(), and v_MultiplyByInvMassMatrix().